SuperNOVAS C++ API v1.6
High-precision C/C++ astrometry library
Loading...
Searching...
No Matches
supernovas::Equatorial Class Reference

Equatorial coordinates (RA, Dec = α, δ), representing the direction ob the sky, for a particular type of equatorial coordinate reference system, relative to the equator and equinox on that system. More...

#include <supernovas.h>

Inheritance diagram for supernovas::Equatorial:

Public Member Functions

 Equatorial (const Angle &ra, const Angle &dec, const Equinox &system=Equinox::icrs())
 Instantiates equatorial coordinates with the specified right-ascention (R.A.) and declination coordinates, optionally specifying a system and a distance if needed.
 Equatorial (const Position &pos, const Equinox &system=Equinox::icrs())
 Instantiates equatorial coordinates with the specified rectangular components.
 Equatorial (const std::string &ra, const std::string &dec, const Equinox &system=Equinox::icrs())
 Instantiates equatorial coordinates with the specified string representations of right-ascention (R.A.) and declination, optionally specifying a system and a distance if needed.
 Equatorial (double ra_rad, double dec_rad, const Equinox &system=Equinox::icrs())
 Instantiates equatorial coordinates with the specified right-ascention (R.A.) and declination coordinates, optionally specifying a system and a distance if needed.
const Angledec () const
 Returns the declination coordinate as an angle.
Angle distance_to (const Equatorial &other) const
 Returns the angular distance of these equatorial coordiantes to/from the specified other equatorial coordinates.
bool equals (const Equatorial &other, const Angle &precision) const
 Checks if these equatorial coordinates are the same as another, within the specified precision.
bool equals (const Equatorial &other, double precision_rad=1.0 *Unit::uas) const
 Checks if these equatorial coordinates are the same as another, within the specified precision.
bool operator!= (const Equatorial &other) const
 Checks if these equatorial coordinates differ from another, by more than 1 μas.
bool operator== (const Equatorial &other) const
 Checks if these equatorial coordinates are the same as another, within 1 μas.
Equatorial operator>> (const Equinox &system) const
 Converts these equatorial coordinates to another equatorial coordinate system.
TimeAngle ra () const
 Returns the right ascention (R.A.) coordinate as a time-angle.
const Equinoxsystem () const
 Returns the equatorial system (type and epoch) in which these equatorial coordinates are defined.
enum novas_reference_system system_type () const
 Retuens the equatorial reference system type in which thse equatorial coordinates are defined.
Equatorial to_cirs (const Time &time) const
 Converts these equatorial coordinates to the Celestial Intermediate Reference System (CIRS) coordinate system, at the specified coordinate epoch.
Ecliptic to_ecliptic () const
 Returns the equivalent ecliptic coordinates corresponding to these equatorial coordinates.
Galactic to_galactic () const
 Returns the equivalent galactic coordinates corresponding to these equatorial coordinates.
Equatorial to_hip () const
 Converts these equatorial coordinates to the Hipparcos catalog coordinate system (= J1991.25).
Equatorial to_icrs () const
 Converts these equatorial coordinates to the International Celestial Reference System (ICRS).
Equatorial to_j2000 () const
 Converts these equatorial coordinates to the J2000 (= FK5) catalog coordinate system.
Equatorial to_mod (const Time &time) const
 Converts these equatorial coordinates to the Mean-of-Date (MOD) catalog coordinate system, at the specified coordinate epoch.
Equatorial to_mod_at_besselian_epoch (double year) const
 Converts these equatorial coordinates to the Mean-of-Date (MOD) catalog coordinate system, at the specified Besselian coordinate epoch.
std::string to_string (enum novas_separator_type separator=NOVAS_SEP_UNITS_AND_SPACES, int decimals=3) const override
 Returns a string representation of these equatorial coordinates in HMS / DMS format, optionally specifying the type of separator to use and the precision to print.
Equatorial to_system (const Equinox &system) const
 Converts these equatorial coordinates to another equatorial coordinate system.
Equatorial to_tod (const Time &time) const
 Converts these equatorial coordinates to the True-of-Date (TOD) coordinate system, at the specified coordinate epoch.
Public Member Functions inherited from supernovas::Spherical
 Spherical (const Angle &longitude, const Angle &latitude)
 Instantiates new spherical coordinates with the specified components.
 Spherical (const std::string &longitude, const std::string &latitude)
 Instantiates spherical coordinates with the specified string representations of the longitude and latitude coordinates, optionally specifying a system and a distance if needed.
 Spherical (double longitude_rad, double latitude_rad)
 Instantiates new spherical coordinates with the specified components.
virtual ~Spherical ()
const Anglelatitude () const
 Returns the latitude coordinate as an angle.
const Anglelongitude () const
 Returns the longitude coordinate as an angle.
Position xyz (const Coordinate &distance) const
 Returns the cartesian position vector corresponding to these spherical coordinates.
Public Member Functions inherited from supernovas::Validating
bool is_valid () const
 Returns the previously set 'valid' stae of the implementing instance.
 operator bool () const
 Objects that implement Validating can be used in conditionals directly, without explicitly calling is_valid().

Static Public Member Functions

static const Equatorial & undefined ()
 Returns a reference to a statically defined standard invalid equatorial coordinates.

Additional Inherited Members

Protected Member Functions inherited from supernovas::Spherical
 Spherical ()
 Instantiates invalid spherical coordinates.
Angle distance_to (const Spherical &other) const
 Returns the angular distance of these spherical coordiantes to/from the specified other spherical coordinates.
bool equals (const Spherical &other, double precision) const
Protected Member Functions inherited from supernovas::Validating
 Validating ()
 dummy constructor;
Protected Attributes inherited from supernovas::Validating
bool _valid = false
 the state variable.

Detailed Description

Equatorial coordinates (RA, Dec = α, δ), representing the direction ob the sky, for a particular type of equatorial coordinate reference system, relative to the equator and equinox on that system.

Constructor & Destructor Documentation

◆ Equatorial() [1/4]

supernovas::Equatorial::Equatorial ( double ra_rad,
double dec_rad,
const Equinox & system = Equinox::icrs() )

Instantiates equatorial coordinates with the specified right-ascention (R.A.) and declination coordinates, optionally specifying a system and a distance if needed.

Parameters
ra_rad[rad] right ascention (R.A.) coordinate
dec_rad[rad] declination coordinate
system(optional) the equatorial coordinate reference system in which the coordinates are specified (default: ICRS)

References supernovas::Spherical::Spherical(), and system().

◆ Equatorial() [2/4]

supernovas::Equatorial::Equatorial ( const Angle & ra,
const Angle & dec,
const Equinox & system = Equinox::icrs() )

Instantiates equatorial coordinates with the specified right-ascention (R.A.) and declination coordinates, optionally specifying a system and a distance if needed.

Parameters
raright ascention (R.A.) coordinate
decdeclination coordinate
system(optional) the equatorial coordinate reference system in which the coordinates are specified (default: ICRS)

References supernovas::Spherical::Spherical(), dec(), ra(), and system().

◆ Equatorial() [3/4]

supernovas::Equatorial::Equatorial ( const std::string & ra,
const std::string & dec,
const Equinox & system = Equinox::icrs() )

Instantiates equatorial coordinates with the specified string representations of right-ascention (R.A.) and declination, optionally specifying a system and a distance if needed.

After instantiation, you should check that the resulting coordinates are valid, e.g. as:

Equatorial coords = Equatorial(..., ...);
if(!coords.is_valid()) {
// oops, looks like the angles could not be parsed...
return;
}
bool is_valid() const
Returns the previously set 'valid' stae of the implementing instance.
Definition supernovas.h:247
Parameters
rastring representation of the right ascention (R.A.) coordinate in HMS or decimal hours.
decstring representation of the declination coordinate as DMS or decimal degrees.
system(optional) the equatorial coordinate reference system in which the coordinates are specified (default: ICRS)
See also
novas_str_hours(), novas_str_degrees() for details on string representation that can be parsed.
novas_parse_hours(), novas_parse_degrees() for more managed parsing from strings.

References dec(), ra(), and system().

◆ Equatorial() [4/4]

supernovas::Equatorial::Equatorial ( const Position & pos,
const Equinox & system = Equinox::icrs() )
explicit

Instantiates equatorial coordinates with the specified rectangular components.

Parameters
posposition vector
system(optional) the equatorial coordinate reference system in which the coordinates are specified (default: ICRS)

References supernovas::Spherical::Spherical(), and system().

Member Function Documentation

◆ dec()

const Angle & supernovas::Equatorial::dec ( ) const

Returns the declination coordinate as an angle.

Returns
the declination coordinate.
See also
ra()

References supernovas::Spherical::latitude().

Referenced by supernovas::CatalogEntry::CatalogEntry(), Equatorial(), Equatorial(), to_ecliptic(), to_galactic(), to_string(), and to_system().

◆ distance_to()

Angle supernovas::Equatorial::distance_to ( const Equatorial & other) const

Returns the angular distance of these equatorial coordiantes to/from the specified other equatorial coordinates.

Parameters
otherthe reference equatorial coordinates
Returns
the angular distance of thereturn Angle::operator+(r);se coordinates to/from the argument.

References supernovas::Spherical::distance_to(), and supernovas::Validating::is_valid().

◆ equals() [1/2]

bool supernovas::Equatorial::equals ( const Equatorial & other,
const Angle & precision ) const

Checks if these equatorial coordinates are the same as another, within the specified precision.

Parameters
otherthe reference equatorial coordinates
precision(optional) precision for equality test (default: 1 μas).
Returns
true if these coordinates are the same as the reference within the precision, or else false.
See also
operator==()

References equals(), and supernovas::Angle::rad().

◆ equals() [2/2]

bool supernovas::Equatorial::equals ( const Equatorial & other,
double precision_rad = 1.0 * Unit::uas ) const

Checks if these equatorial coordinates are the same as another, within the specified precision.

Parameters
otherthe reference equatorial coordinates
precision_rad[rad] (optional) precision for equality test (default: 1 μas).
Returns
true if these coordinates are the same as the reference within the precision, or else false.
See also
equals(), operator==()

References supernovas::Spherical::equals().

Referenced by equals(), operator!=(), and operator==().

◆ operator!=()

bool supernovas::Equatorial::operator!= ( const Equatorial & other) const

Checks if these equatorial coordinates differ from another, by more than 1 μas.

Parameters
otherthe reference equatorial coordinates
Returns
true if these coordinates differ from the reference, by more than 1 μas, or else false.
See also
operator==()

References equals().

◆ operator==()

bool supernovas::Equatorial::operator== ( const Equatorial & other) const

Checks if these equatorial coordinates are the same as another, within 1 μas.

Parameters
otherthe reference equatorial coordinates
Returns
true if these coordinates are the same as the reference within 1 μas, or else false.
See also
operator!=()

References equals().

◆ operator>>()

Equatorial supernovas::Equatorial::operator>> ( const Equinox & system) const

Converts these equatorial coordinates to another equatorial coordinate system.

Same as to_system().

Parameters
systemthe equatorial coordinate system (type and epoch) to convert to.
Returns
new equatorial coordinates, which represent the same equatorial position as this, but expressed in the specified other coordinate reference system.
See also
to_system()

References supernovas::Validating::is_valid(), system(), and to_system().

◆ ra()

TimeAngle supernovas::Equatorial::ra ( ) const

Returns the right ascention (R.A.) coordinate as a time-angle.

Returns
the right ascention (R.A.) coordinate.
See also
dec()

References supernovas::Validating::is_valid(), and supernovas::Spherical::longitude().

Referenced by supernovas::CatalogEntry::CatalogEntry(), Equatorial(), Equatorial(), to_ecliptic(), to_galactic(), to_string(), and to_system().

◆ system()

const Equinox & supernovas::Equatorial::system ( ) const

Returns the equatorial system (type and epoch) in which these equatorial coordinates are defined.

Returns
the coordinate reference system (type and epoch).
See also
system_type()

Referenced by Equatorial(), Equatorial(), Equatorial(), Equatorial(), operator>>(), and to_system().

◆ system_type()

enum novas_reference_system supernovas::Equatorial::system_type ( ) const

Retuens the equatorial reference system type in which thse equatorial coordinates are defined.

Returns
the type of coordinate reference system
See also
system()

◆ to_cirs()

Equatorial supernovas::Equatorial::to_cirs ( const Time & time) const

Converts these equatorial coordinates to the Celestial Intermediate Reference System (CIRS) coordinate system, at the specified coordinate epoch.

CIRS is defined on the true dynamical equator of date, with its origin at the Celestial Intermediate Origin (CIO).

Parameters
Time of observation and astronomical timescales[day] the astronomical time specification for the coordinate epoch.
Returns
new equatorial coordinates, which represent the same equatorial position as this, but with respect to the true equator and CIO of date.
See also
to_system(), to_tod(), to_icrs()

References supernovas::Equinox::cirs(), supernovas::Validating::is_valid(), and to_system().

◆ to_hip()

Equatorial supernovas::Equatorial::to_hip ( ) const

Converts these equatorial coordinates to the Hipparcos catalog coordinate system (= J1991.25).

Returns
new equatorial coordinates, which represent the same equatorial position as this, but expressed in the Hipparcos (= J1991.25) catalog system.
See also
to_system(), to_icrs(), to_j2000()

References supernovas::Equinox::hip(), supernovas::Validating::is_valid(), and to_system().

◆ to_j2000()

Equatorial supernovas::Equatorial::to_j2000 ( ) const

Converts these equatorial coordinates to the J2000 (= FK5) catalog coordinate system.

Returns
new equatorial coordinates, which represent the same equatorial position as this, but expressed in the J2000 (= FK5) catalog system.
See also
to_system(), to_icrs(), to_hip(), to_mod(), to_mod_at_besselian_epoch(), to_tod(), to_cirs()

References supernovas::Validating::is_valid(), supernovas::Equinox::j2000(), and to_system().

Referenced by supernovas::Ecliptic::to_j2000().

◆ to_mod()

Equatorial supernovas::Equatorial::to_mod ( const Time & time) const

Converts these equatorial coordinates to the Mean-of-Date (MOD) catalog coordinate system, at the specified coordinate epoch.

Parameters
Time of observation and astronomical timescales[day] the astronomical time specification for the coordinate epoch.
Returns
new equatorial coordinates, which represent the same equatorial position as this, but expressed in the MOD catalog system of date.
See also
to_mod_at_besselian_epoch(), to_system(), to_j2000(), to_tod()

References supernovas::Validating::is_valid(), supernovas::Equinox::mod(), and to_system().

Referenced by supernovas::Ecliptic::to_mod().

◆ to_mod_at_besselian_epoch()

Equatorial supernovas::Equatorial::to_mod_at_besselian_epoch ( double year) const

Converts these equatorial coordinates to the Mean-of-Date (MOD) catalog coordinate system, at the specified Besselian coordinate epoch.

Parameters
year[yr] Besselian year for the coordinate epoch (e.g. 1950.0 for B1950).
Returns
new equatorial coordinates, which represent the same equatorial position as this, but expressed in the catalog system of the specified Besselian epoch.
See also
to_mod(), to_system(), to_j2000(), to_tod()

References supernovas::Validating::is_valid(), supernovas::Equinox::mod_at_besselian_epoch(), and to_system().

◆ to_string()

std::string supernovas::Equatorial::to_string ( enum novas_separator_type separator = NOVAS_SEP_UNITS_AND_SPACES,
int decimals = 3 ) const
overridevirtual

Returns a string representation of these equatorial coordinates in HMS / DMS format, optionally specifying the type of separator to use and the precision to print.

Parameters
separator(optional) the type of separators to use for the HMS / DMS representations of the components (default: units and spaces).
decimals(optional) the number of decimal places to print for the seconds (default: 3)
Returns
a new string with the human-readable representation of these equatorial coordinates.

Reimplemented from supernovas::Spherical.

References dec(), ra(), supernovas::Angle::to_string(), and supernovas::TimeAngle::to_string().

Referenced by supernovas::CatalogEntry::to_string().

◆ to_tod()

Equatorial supernovas::Equatorial::to_tod ( const Time & time) const

Converts these equatorial coordinates to the True-of-Date (TOD) coordinate system, at the specified coordinate epoch.

TOD is defined on the true dynamical equator of date, with its origin at the true equinox of date.

Parameters
Time of observation and astronomical timescales[day] the astronomical time specification for the coordinate epoch.
Returns
new equatorial coordinates, which represent the same equatorial position as this, but expressed with respect to the true equator and equinox of date.
See also
to_system(), to_cirs(), to_j2000(), to_mod()

References supernovas::Validating::is_valid(), to_system(), and supernovas::Equinox::tod().

Referenced by supernovas::Ecliptic::to_tod().

◆ undefined()

const Equatorial & supernovas::Equatorial::undefined ( )
static

Returns a reference to a statically defined standard invalid equatorial coordinates.

These invalid coordinates may be used inside any object that is invalid itself.

Returns
a reference to a static standard invalid equatorial coordinates.

Referenced by supernovas::EquatorialTrack::projected_at(), supernovas::Ecliptic::to_equatorial(), supernovas::Galactic::to_equatorial(), and to_system().