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

Ecliptic coordinates (l, b or λ, β), representing the direction on the sky, for a particular type of equatorial coordinate reference system, relative to the ecliptic and equinox of that system. More...

#include <supernovas.h>

Inheritance diagram for supernovas::Ecliptic:

Public Member Functions

 Ecliptic (const Angle &longitude, const Angle &latitude, const Equinox &system=Equinox::icrs())
 Instantiates new ecliptic coordinates with the specified parameters.
 Ecliptic (const Position &pos, const Equinox &system=Equinox::icrs())
 Instantiates new ecliptic coordinates with the specified ecliptic cartesian position vector.
 Ecliptic (const std::string &longitude, const std::string &latitude, const Equinox &system=Equinox::icrs())
 Instantiates ecliptic coordinates with the specified string representations of the longitude and latitude coordinates, optionally specifying a system and a distance if needed.
 Ecliptic (double longitude_rad, double latitude_rad, const Equinox &system=Equinox::icrs())
 Instantiates new ecliptic coordinates with the specified parameters.
Angle distance_to (const Ecliptic &other) const
 Returns the angular distance of these ecliptic coordiantes to/from the specified other ecliptic coordinates.
bool equals (const Ecliptic &other, const Angle &precision) const
 Checks if these ecliptic coordinates are the same as another, within the specified precision.
bool equals (const Ecliptic &other, double precision_rad=1.0 *Unit::uas) const
 Checks if these ecliptic coordinates are the same as another, within the specified precision.
enum novas_equator_type equator_type () const
 Returns the type of equator (ICRS, mean, or true) that is used for these ecliptic coordinates.
double jd () const
 Returns the Julian date of the epoch for which the coordinates are defined.
double mjd () const
 Returns the Modified Julian Date (MJD) of the epoch for which the coordinates are defined.
bool operator!= (const Ecliptic &other) const
 Checks if these ecliptic coordinates differ from another, by more than 1 μ uas.
bool operator== (const Ecliptic &other) const
 Checks if these ecliptic coordinates are the same as another, within 1 μas.
Ecliptic operator>> (const Equinox &system) const
 Converts these ecliptic coordinates to the ecliptic coordinate system with respect to the specified equinox of date.
Equinox system () const
 Returns the equinox of date relative to which these ecliptic coordinates are defined.
Equatorial to_equatorial () const
 Returns a reference to a statically defined standard invalid Galactic coordinates.
Galactic to_galactic () const
 Converts these ecliptic coordinates to equivalent Galactic coordinates.
Ecliptic to_icrs () const
 Converts these ecliptic coordinates to ICRS ecliptic coordinates.
Ecliptic to_j2000 () const
 Converts these ecliptic coordinates to J2000 ecliptic coordinates.
Ecliptic to_mod (const Time &time) const
 Converts these ecliptic coordinates to Mean-of-Date (MOD) ecliptic coordinates at the specified epch.
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 ecliptic coordinates, optionally specifying a type of separator to use for the DMS angles, and the decimal places to show for the seconds.
Ecliptic to_system (const Equinox &system) const
 Converts these ecliptic coordinates to the ecliptic coordinate system with respect to the specified equinox of date.
Ecliptic to_tod (const Time &time) const
 Converts these ecliptic coordinates to True-of-Date (TOD) ecliptic coordinates at the specified epch.
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 Ecliptic & undefined ()
 Returns a reference to a statically defined standard invalid ecliptic 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

Ecliptic coordinates (l, b or λ, β), representing the direction on the sky, for a particular type of equatorial coordinate reference system, relative to the ecliptic and equinox of that system.

Constructor & Destructor Documentation

◆ Ecliptic() [1/4]

supernovas::Ecliptic::Ecliptic ( double longitude_rad,
double latitude_rad,
const Equinox & system = Equinox::icrs() )

Instantiates new ecliptic coordinates with the specified parameters.

Parameters
longitude_rad[rad] ecliptic longitude coordinate
latitude_rad[rad] ecliptic latitude coordinate
system(optional) the equatorial coordinate reference system that defines the the origin of ecliptic longitude, that is the equinox of date (default: ICRS).

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

◆ Ecliptic() [2/4]

supernovas::Ecliptic::Ecliptic ( const Angle & longitude,
const Angle & latitude,
const Equinox & system = Equinox::icrs() )

Instantiates new ecliptic coordinates with the specified parameters.

Parameters
longitudeecliptic longitude coor dinate
latitudeecliptic latitude coordinate
system(optional) The equatorial coordinate reference system that defines the origin of ecliptic longitude, that is the equinox of date (default: ICRS).

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

◆ Ecliptic() [3/4]

supernovas::Ecliptic::Ecliptic ( const std::string & longitude,
const std::string & latitude,
const Equinox & system = Equinox::icrs() )

Instantiates ecliptic coordinates with the specified string representations of the longitude and latitude coordinates, optionally specifying a system and a distance if needed.

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

Ecliptic coords = Ecliptic(..., ...);
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
longitudestring representation of the longitude coordinate in DMS or a decimnal degrees.
latitudestring representation of the declination coordinate as DMS or decimal degrees.
system(optional) the equatorial coordinate reference system that defines the the origin of ecliptic longitude, that is the equinox of date (default: ICRS).
See also
novas_str_degrees() for details on string representation that can be parsed.
novas_parse_degrees() for more managed parsing from strings.

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

◆ Ecliptic() [4/4]

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

Instantiates new ecliptic coordinates with the specified ecliptic cartesian position vector.

Parameters
posEcliptic xyz position vector
system(optional) the equatorial coordinate reference system that defines the the origin of ecliptic longitude, that is the equinox of date (default: ICRS).

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

Member Function Documentation

◆ distance_to()

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

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

Parameters
otherthe reference ecliptic coordinates
Returns
the angular distance of these coordinates to/from the argument.

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

◆ equals() [1/2]

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

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

Parameters
otherthe reference ecliptic 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::Ecliptic::equals ( const Ecliptic & other,
double precision_rad = 1.0 * Unit::uas ) const

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

Parameters
otherthe reference ecliptic 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
operator==()

References supernovas::Spherical::equals().

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

◆ equator_type()

enum novas_equator_type supernovas::Ecliptic::equator_type ( ) const

Returns the type of equator (ICRS, mean, or true) that is used for these ecliptic coordinates.

Returns
the type of equator that defines the origin (equinox), uch as ICRS, mean, or true.
See also
Equinox::equator_type()

Referenced by Ecliptic(), and Ecliptic().

◆ jd()

double supernovas::Ecliptic::jd ( ) const

Returns the Julian date of the epoch for which the coordinates are defined.

Returns
[day] the (TDB-based) Julian date of the epoch for which the coordinates are defined.
See also
mjd()

Referenced by Ecliptic(), and Ecliptic().

◆ mjd()

double supernovas::Ecliptic::mjd ( ) const

Returns the Modified Julian Date (MJD) of the epoch for which the coordinates are defined.

Returns
[day] the (TDB-based) MJD of the epoch for which the coordinates are defined.
See also
jd()

References NOVAS_JD_MJD0.

◆ operator!=()

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

Checks if these ecliptic coordinates differ from another, by more than 1 μ uas.

Parameters
otherthe reference ecliptic 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::Ecliptic::operator== ( const Ecliptic & other) const

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

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

References equals().

◆ operator>>()

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

Converts these ecliptic coordinates to the ecliptic coordinate system with respect to the specified equinox of date.

Same as to_system().

Parameters
systemthe requested equinox of date for returned coordinates.
Returns
new ecliptic coordinates, which represent the same position as this, but expressed relaive to the specified equinox.
See also
to_system()

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

◆ system()

Equinox supernovas::Ecliptic::system ( ) const

Returns the equinox of date relative to which these ecliptic coordinates are defined.

Returns
the equinox of date, which these ecliptic coordinates are based on.

References supernovas::Equinox::icrs(), supernovas::Equinox::j2000(), supernovas::Equinox::mod(), NOVAS_GCRS_EQUATOR, NOVAS_JD_J2000, NOVAS_MEAN_EQUATOR, NOVAS_TRUE_EQUATOR, supernovas::Equinox::tod(), and supernovas::Equinox::undefined().

Referenced by Ecliptic(), Ecliptic(), Ecliptic(), Ecliptic(), distance_to(), operator>>(), to_equatorial(), and to_system().

◆ to_icrs()

Ecliptic supernovas::Ecliptic::to_icrs ( ) const

Converts these ecliptic coordinates to ICRS ecliptic coordinates.

Returns
the equivalent ICRS ecliptic coordinates.
See also
to_system(), to_j2000(), to_mod(), to_tod()

References supernovas::Validating::is_valid(), NOVAS_GCRS_EQUATOR, supernovas::Equatorial::to_ecliptic(), to_equatorial(), and supernovas::Equatorial::to_icrs().

◆ to_j2000()

Ecliptic supernovas::Ecliptic::to_j2000 ( ) const

Converts these ecliptic coordinates to J2000 ecliptic coordinates.

Returns
the equivalent J2000 ecliptic coordinates.
See also
to_system(), to_icrs(), to_mod(), to_tod()

References supernovas::Validating::is_valid(), NOVAS_JD_J2000, NOVAS_MEAN_EQUATOR, supernovas::Equatorial::to_ecliptic(), to_equatorial(), and supernovas::Equatorial::to_j2000().

Referenced by to_mod().

◆ to_mod()

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

Converts these ecliptic coordinates to Mean-of-Date (MOD) ecliptic coordinates at the specified epch.

Parameters
Time of observation and astronomical timescalesthe astronomical time specifying the coordinate epoch.
Returns
the equivalent MOD ecliptic coordinates at the specified date.
See also
to_system(), to_mod(), to_tod(), to_icrs(), to_j2000()

References supernovas::Validating::is_valid(), supernovas::Time::j2000(), supernovas::Time::jd(), NOVAS_MEAN_EQUATOR, supernovas::Equatorial::to_ecliptic(), to_equatorial(), to_j2000(), and supernovas::Equatorial::to_mod().

◆ to_string()

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

Returns a string representation of these ecliptic coordinates, optionally specifying a type of separator to use for the DMS angles, and the decimal places to show for the seconds.

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

Reimplemented from supernovas::Spherical.

References supernovas::Spherical::to_string().

◆ to_system()

Ecliptic supernovas::Ecliptic::to_system ( const Equinox & system) const

Converts these ecliptic coordinates to the ecliptic coordinate system with respect to the specified equinox of date.

Parameters
systemthe requested equinox of date for returned coordinates.
Returns
new ecliptic coordinates, which represent the same position as this, but expressed relaive to the specified equinox.
See also
operator>>(), to_icrs(), to_j2000(), to_mod(), to_tod()

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

Referenced by operator>>().

◆ to_tod()

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

Converts these ecliptic coordinates to True-of-Date (TOD) ecliptic coordinates at the specified epch.

Parameters
Time of observation and astronomical timescalesthe astronomical time specifying the coordinate epoch.
Returns
the equivalent TOD ecliptic coordinates at the specified date.
See also
to_system(), to_tod(), to_mod(), to_icrs(), to_j2000()

References supernovas::Validating::is_valid(), supernovas::Time::jd(), NOVAS_TRUE_EQUATOR, supernovas::Equatorial::to_ecliptic(), to_equatorial(), and supernovas::Equatorial::to_tod().

◆ undefined()

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

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

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

Returns
a reference to the static standard invalid coordinates.

Referenced by supernovas::Equatorial::to_ecliptic().