CppEphem
|
CECoordinates class is responsible for doing all of the celestial coordinate conversions. Command line tools are also available for "on the fly" conversions.
The different coordinate systems used:
Definition at line 48 of file CECoordinates.h.
Public Member Functions | |
virtual CEAngle | AngularSeparation (const CECoordinates &coords) const |
Get the angular separation between the coordinates represented by this object and another coordinate object. More... | |
CECoordinates () | |
Default constructor. More... | |
CECoordinates (const CEAngle &xcoord, const CEAngle &ycoord, const CECoordinateType &coord_type=CECoordinateType::ICRS) | |
Primary constructor (NOTE: xcoord & ycoord are expected to be in radians by default. More... | |
CECoordinates (const CECoordinates &other) | |
Copy constructor. More... | |
CECoordinates (const CECoordinateType &coord_type) | |
Constructor from a coordinate type. More... | |
CECoordinates (const std::vector< double > &xcoord, const std::vector< double > &ycoord, const CECoordinateType &coord_type=CECoordinateType::ICRS) | |
Primary constructor. More... | |
CECoordinates | ConvertTo (CECoordinateType output_coord_type, const CEObserver &observer, const CEDate &date=CEDate::CurrentJD()) |
Convert these coordinates to another coordinate system NOTE: If this object is not OBSERVED coordinates, only JD is needed. More... | |
CECoordinates | ConvertTo (CECoordinateType output_coord_type, double jd=CEDate::CurrentJD(), double longitude=0.0, double latitude=0.0, double elevation_m=0.0, double pressure_hPa=-1.0, double temperature_celsius=-1000, double relative_humidity=0.0, double dut1=0.0, double xp=0.0, double yp=0.0, double wavelength_um=0.5) |
Convert these coordinates to another coordinate system NOTE: If neither this object or output_coord_type are OBSERVED coordinates, only JD is needed. More... | |
CECoordinates | ConvertToCIRS (double jd=CEDate::CurrentJD(), double longitude=0.0, double latitude=0.0, double elevation_m=0.0, double pressure_hPa=-1.0, double temperature_celsius=-1000, double relative_humidity=0.0, double dut1=0.0, double xp=0.0, double yp=0.0, double wavelength_um=0.5) |
Convert these coordinates to CIRS coordinates. More... | |
CECoordinates | ConvertToGalactic (double jd=CEDate::CurrentJD(), double longitude=0.0, double latitude=0.0, double elevation_m=0.0, double pressure_hPa=-1.0, double temperature_celsius=-1000, double relative_humidity=0.0, double dut1=0.0, double xp=0.0, double yp=0.0, double wavelength_um=0.5) |
Convert these coordinates to GALACTIC coordinates. More... | |
CECoordinates | ConvertToICRS (double jd=CEDate::CurrentJD(), double longitude=0.0, double latitude=0.0, double elevation_m=0.0, double pressure_hPa=-1.0, double temperature_celsius=-1000, double relative_humidity=0.0, double dut1=0.0, double xp=0.0, double yp=0.0, double wavelength_um=0.5) |
Convert these coordinates to ICRS coordinates. More... | |
CECoordinates | ConvertToObserved (double jd=CEDate::CurrentJD(), double longitude=0.0, double latitude=0.0, double elevation_m=0.0, double pressure_hPa=-1.0, double temperature_celsius=-1000, double relative_humidity=0.0, double dut1=0.0, double xp=0.0, double yp=0.0, double wavelength_um=0.5) |
Convert these coordinates to observed coordinates. More... | |
CECoordinateType | GetCoordSystem (void) const |
Return coordinate system. More... | |
virtual CECoordinates | GetObservedCoords (const CEDate &date, const CEObserver &observer) const |
Return the observed coordinates using an observer object (CEObserver) More... | |
virtual CECoordinates | GetObservedCoords (const double &julian_date, const double &longitude, const double &latitude, const double &elevation_m=0.0, const double &pressure_hPa=-1.0, const double &temperature_celsius=-1000, const double &relative_humidity=0.0, const double &dut1=0.0, const double &xp=0.0, const double &yp=0.0, const double &wavelength_um=0.5) const |
Return the local sky coordinates of this object as a CECoordinates object. More... | |
CECoordinates & | operator= (const CECoordinates &other) |
Overloaded '=' (assignment) operator. More... | |
std::string | print (void) const |
Generate a message string that specifies the information about this coordinate. More... | |
virtual void | SetCoordinates (const CEAngle &xcoord, const CEAngle &ycoord, const CECoordinateType &coord_type=CECoordinateType::ICRS) |
Set the coordinates of this object. More... | |
virtual void | SetCoordinates (const CECoordinates &coords) |
Set the coordinates from another CECoordinates object. More... | |
virtual CEAngle | XCoord (const double &jd=CppEphem::julian_date_J2000()) const |
Return x coordinate at given Julian date. More... | |
virtual double | XCoordinate_Deg (double jd=CppEphem::julian_date_J2000()) const |
Returns x coordinate at given julian date in degrees. More... | |
virtual double | XCoordinate_Rad (double jd=CppEphem::julian_date_J2000()) const |
Return x coordinate at given Julian date in radians. More... | |
virtual CEAngle | YCoord (const double &jd=CppEphem::julian_date_J2000()) const |
Return y coordinate at given Julian date. More... | |
virtual double | YCoordinate_Deg (double jd=CppEphem::julian_date_J2000()) const |
Returns y coordinate at given Julian date in degrees. More... | |
virtual double | YCoordinate_Rad (double jd=CppEphem::julian_date_J2000()) const |
Returns y coordinate at given Julian date in radians. More... | |
virtual | ~CECoordinates () |
Destructor. More... | |
Static Public Member Functions | |
static CEAngle | AngularSeparation (const CEAngle &xcoord_first, const CEAngle &ycoord_first, const CEAngle &xcoord_second, const CEAngle &ycoord_second) |
Get the angular separation between two sets of coordinates. More... | |
static CEAngle | AngularSeparation (const CECoordinates &coords1, const CECoordinates &coords2) |
Get the angular separation between two coordinate objects. More... | |
static void | CIRS2Galactic (double ra, double dec, double *glon, double *glat, const CEDate &date=CEDate(), const CEAngleType &angle_type=CEAngleType::RADIANS) |
CIRS -> Galactic coordinate conversion. More... | |
static void | CIRS2ICRS (double input_ra, double input_dec, double *return_ra, double *return_dec, const CEDate &date=CEDate(), const CEAngleType &angle_type=CEAngleType::RADIANS) |
CIRS -> ICRS coordinate conversion. More... | |
static int | CIRS2Observed (double ra, double dec, double *az, double *zen, const CEDate &date, const CEObserver &observer, const CEAngleType &angle_type=CEAngleType::RADIANS, double wavelength_um=0.5, double *observed_ra=nullptr, double *observed_dec=nullptr, double *hour_angle=nullptr) |
CIRS -> Observed (or observer specific) coordinate conversion This function takes in verious observation parameters The integer returned is a status code with the following meanings: +1 = dubious year (too far into the past/future to be trusted) 0 = OK status -1 = unacceptable date. More... | |
static int | CIRS2Observed (double ra, double dec, double *az, double *zen, double julian_date, double longitude, double latitude, double elevation_m=0.0, double pressure_hPa=-1.0, double temperature_celsius=-1000, double relative_humidity=0.0, double dut1=0.0, double xp=0.0, double yp=0.0, double wavelength_um=0.5, double *observed_ra=nullptr, double *observed_dec=nullptr, double *hour_angle=nullptr) |
Raw method for converting CIRS -> Observed (observer specific) coordinates (uses the SOFA 'iauAtio13' function) Note: All angles are expected to be in radians. More... | |
static double | DMSToAngle (const std::vector< double > &angle, const CEAngleType &return_type=CEAngleType::DEGREES) |
Convert a given angle vector from {degrees, minutes, seconds} to an angle. More... | |
static void | Galactic2CIRS (double glon, double glat, double *ra, double *dec, const CEDate &date=CEDate(), const CEAngleType &angle_type=CEAngleType::RADIANS) |
Galactic -> CIRS coordinate conversion. More... | |
static void | Galactic2ICRS (double glon, double glat, double *ra, double *dec, const CEAngleType &angle_type=CEAngleType::RADIANS) |
Galactic -> ICRS coordinate conversion. More... | |
static int | Galactic2Observed (double glon, double glat, double *az, double *zen, const CEDate &date, const CEObserver &observer, const CEAngleType &angle_type=CEAngleType::RADIANS, double wavelength_um=0.5, double *observed_glon=nullptr, double *observed_glat=nullptr, double *hour_angle=nullptr) |
Galactic -> Observed (i.e. More... | |
static int | Galactic2Observed (double glon, double glat, double *az, double *zen, double julian_date, double longitude, double latitude, double elevation_m=0.0, double pressure_hPa=-1.0, double temperature_celsius=-1000, double relative_humidity=0.0, double dut1=0.0, double xp=0.0, double yp=0.0, double wavelength_um=0.50, double *observed_glon=nullptr, double *observed_glat=nullptr, double *hour_angle=nullptr) |
Raw method for converting Galactic -> Observed (observer specific) coordinates. More... | |
static std::vector< double > | GetDMS (const double &angle, const CEAngleType &angle_type=CEAngleType::DEGREES) |
Convert a given angle into degrees, arcminutes, arcseconds. More... | |
static std::vector< double > | GetHMS (const double &angle, const CEAngleType &angle_type=CEAngleType::DEGREES) |
Convert a given angle into hours, minutes, seconds. More... | |
static double | HMSToAngle (const std::vector< double > &angle, const CEAngleType &return_type=CEAngleType::DEGREES) |
Convert from {hours, minutes, seconds} to an angle. More... | |
static void | ICRS2CIRS (double input_ra, double input_dec, double *return_ra, double *return_dec, const CEDate &date=CEDate(), const CEAngleType &angle_type=CEAngleType::RADIANS) |
ICRS -> CIRS coordinate conversion. More... | |
static void | ICRS2Galactic (double ra, double dec, double *glon, double *glat, const CEAngleType &angle_type=CEAngleType::RADIANS) |
ICRS -> Galactic coordinate conversion (uses the SOFA 'iauIcrs2g' function) More... | |
static int | ICRS2Observed (double ra, double dec, double *az, double *zen, const CEDate &date, const CEObserver &observer, const CEAngleType &angle_type=CEAngleType::RADIANS, double wavelength_um=0.5, double *observed_ra=nullptr, double *observed_dec=nullptr, double *hour_angle=nullptr) |
ICRS -> Observed coordinate conversion. More... | |
static int | ICRS2Observed (double ra, double dec, double *az, double *zen, double julian_date, double longitude, double latitude, double elevation_m=0.0, double pressure_hPa=-1.0, double temperature_celsius=-1000, double relative_humidity=0.0, double dut1=0.0, double xp=0.0, double yp=0.0, double wavelength_um=0.5, double *observed_ra=nullptr, double *observed_dec=nullptr, double *hour_angle=nullptr) |
Raw method for converting CIRS -> Observed (observer specific) coordinates. More... | |
static int | Observed2CIRS (double az, double zen, double *ra, double *dec, const CEDate &date, const CEObserver &observer, const CEAngleType &angle_type=CEAngleType::RADIANS) |
Convert Observed -> CIRS coordinates. More... | |
static int | Observed2CIRS (double az, double zen, double *ra, double *dec, double julian_date, double longitude, double latitude, double elevation_m=0.0, double pressure_hPa=-1.0, double temperature_celsius=-1000, double relative_humidity=0.0, double dut1=0.0, double xp=0.0, double yp=0.0, double wavelength_um=0.5) |
Raw method for converting Observed (observer specific) -> CIRS coordinates (uses the SOFA 'iauAtoi13' function) Note: All angles are expected to be in radians. More... | |
static int | Observed2Galactic (double az, double zen, double *glon, double *glat, const CEDate &date, const CEObserver &observer, const CEAngleType &angle_type=CEAngleType::RADIANS) |
Convert Observed -> Galactic coordinates. More... | |
static int | Observed2Galactic (double az, double zen, double *glon, double *glat, double julian_date, double longitude, double latitude, double elevation_m=0.0, double pressure_hPa=-1.0, double temperature_celsius=-1000, double relative_humidity=0.0, double dut1=0.0, double xp=0.0, double yp=0.0, double wavelength_um=0.50) |
Raw method for converting Observed (observer specific) -> Galactic coordinates Note: All angles are expected to be in radians. More... | |
static int | Observed2ICRS (double az, double zen, double *ra, double *dec, const CEDate &date, const CEObserver &observer, const CEAngleType &angle_type=CEAngleType::RADIANS) |
Convert Observed -> ICRS coordinates. More... | |
static int | Observed2ICRS (double az, double zen, double *ra, double *dec, double julian_date, double longitude, double latitude, double elevation_m=0.0, double pressure_hPa=-1.0, double temperature_celsius=-1000, double relative_humidity=0.0, double dut1=0.0, double xp=0.0, double yp=0.0, double wavelength_um=0.5) |
Raw method for converting Observed (observer specific) -> ICRS coordinates Note: All angles are expected to be in radians. More... | |
Protected Attributes | |
CECoordinateType | coord_type_ |
CEAngle | xcoord_ |
CEAngle | ycoord_ |
Private Member Functions | |
void | copy_members (const CECoordinates &other) |
Copy data members from another CECoordinates object. More... | |
void | free_members (void) |
Cleanup data members that need to be freed or cleared. More... | |
void | init_members (void) |
Set initial values and allocate memory for data members. More... | |
Friends | |
bool | operator!= (const CECoordinates &lhs, const CECoordinates &rhs) |
Compare two coordinate objects. More... | |
bool | operator== (const CECoordinates &lhs, const CECoordinates &rhs) |
Compare two coordinate objects. More... | |
#include <CECoordinates.h>
CECoordinates::CECoordinates | ( | ) |
Default constructor.
Definition at line 38 of file CECoordinates.cpp.
CECoordinates::CECoordinates | ( | const CEAngle & | xcoord, |
const CEAngle & | ycoord, | ||
const CECoordinateType & | coord_type = CECoordinateType::ICRS |
||
) |
Primary constructor (NOTE: xcoord & ycoord are expected to be in radians by default.
[in] | xcoord | X-Coordinate (radians) |
[in] | ycoord | Y-Coordinate (radians) |
[in] | coord_type | Coordinate type (see CECoordinateType) |
[in] | angle_type | Angle type (either DEGREES or RADIANS) |
Definition at line 52 of file CECoordinates.cpp.
CECoordinates::CECoordinates | ( | const std::vector< double > & | xcoord, |
const std::vector< double > & | ycoord, | ||
const CECoordinateType & | coord_type = CECoordinateType::ICRS |
||
) |
Primary constructor.
[in] | xcoord | X-Coordinate {hours,minutes,seconds} |
[in] | ycoord | Y-Coordinate {degrees, arcmin, arcsec} |
[in] | coord_type | Coordinate type (see CECoordinateType) |
[in] | angle_type | Angle type (either DEGREES or RADIANS) |
Definition at line 69 of file CECoordinates.cpp.
CECoordinates::CECoordinates | ( | const CECoordinateType & | coord_type | ) |
Constructor from a coordinate type.
[in] | coord_type | Coordinate type (see CECoordinateType) |
Definition at line 91 of file CECoordinates.cpp.
CECoordinates::CECoordinates | ( | const CECoordinates & | other | ) |
Copy constructor.
[in] | other | Coordinates object to be copied |
Definition at line 102 of file CECoordinates.cpp.
|
virtual |
Destructor.
Definition at line 112 of file CECoordinates.cpp.
|
static |
Get the angular separation between two sets of coordinates.
Method utilizes the SOFA 'iauSeps' algorithm. NOTE: The coordinates are both expected to be in the same coordinate system!
[in] | xcoord_first | X-value for first set of coordinates |
[in] | ycoord_first | Y-value for first set of coordinates |
[in] | xcoord_second | X-value for second set of coordiantes |
[in] | ycoord_second | Y-value for second set of coordinates |
[in] | return_angle_type | Specify whether input angles are DEGREES or RADIANS. (output angle will be in the same format) |
Note that the x-coordinates are expected in the range [0, 2pi] and the y-coordinates are expected in the range [-pi, pi]. Because of this, you need to pass:
Definition at line 1235 of file CECoordinates.cpp.
|
virtual |
Get the angular separation between the coordinates represented by this object and another coordinate object.
NOTE: The coordinates are both expected to be in the same coordinate system! If they are in different coordinate systems, use "ConvertTo()" first.
[in] | coords | Another set of coordinates |
Definition at line 1169 of file CECoordinates.cpp.
|
static |
Get the angular separation between two coordinate objects.
NOTE: The coordinates are both expected to be in the same coordinate system! If they are in different coordinate systems, use "ConvertTo()" first.
[in] | coords1 | First set of coordinates |
[in] | coords2 | Second set of coordinates |
[in] | return_angle_type | Specify whether to return angle as DEGREES or RADIANS |
Note that the x-coordinates are expected in the range [0, 2pi] and the y-coordinates are expected in the range [-pi, pi]. Because of this, OBSERVED coordinates first convert the zenith angle to altitude
Definition at line 1189 of file CECoordinates.cpp.
|
static |
CIRS -> Galactic coordinate conversion.
[in] | input_ra | CIRS right ascension |
[in] | input_dec | CIRS declination |
[out] | glon | Galactic longitude |
[out] | glat | Galactic latitude |
[in] | date | Date information |
[in] | angle_type | Angle format (DEGREES or RADIANS) |
Definition at line 181 of file CECoordinates.cpp.
|
static |
CIRS -> ICRS coordinate conversion.
(uses SOFA 'iauAtic13' function)
[in] | input_ra | CIRS right ascension |
[in] | input_dec | CIRS declination |
[out] | return_ra | ICRS right ascension (returned) |
[out] | return_dec | ICRS declinaton (returned) |
[in] | date | Date information |
[in] | angle_type | Input/output angle type (DEGREES or RADIANS) |
Definition at line 144 of file CECoordinates.cpp.
|
static |
CIRS -> Observed (or observer specific) coordinate conversion This function takes in verious observation parameters The integer returned is a status code with the following meanings: +1 = dubious year (too far into the past/future to be trusted) 0 = OK status -1 = unacceptable date.
[in] | ra | Right ascension in CIRS coordinates |
[in] | dec | Declination in CIRS coordinates |
[out] | az | Azimuth (returned) |
[out] | zen | Zenith angle (returned) |
[in] | observer | Observer information |
[in] | angle_type | Angle type (see CEAngleType) |
[in] | wavelength | Wavelength of the light being observed (micrometers) |
[out] | observed_ra | Observed right ascension (returned) |
[out] | observed_dec | Observed declination (returned) |
[out] | hour_angle | Hour angle of object being observed (returned) |
Definition at line 224 of file CECoordinates.cpp.
|
static |
Raw method for converting CIRS -> Observed (observer specific) coordinates (uses the SOFA 'iauAtio13' function) Note: All angles are expected to be in radians.
[in] | ra | CIRS right ascension (radians) |
[in] | dec | CIRS declination (radians) |
[out] | az | Observed azimuth angle (radians, returned) |
[out] | zen | Observed zenith angle (radians, returned) |
[in] | julian_date | Julian date for conversion |
[in] | longitude | Observer geographic longitude (radians) |
[in] | latitude | Observer geographic latitude (radians) |
[in] | elevation_m | Observer elevation (meters) |
[in] | pressure_hPa | Atmospheric pressure (HPa) |
[in] | temperature_celsius | Temperature (degrees Celsius) |
[in] | relative_humidity | Relative humidity (0.0 - 1.0) |
[in] | dut1 | UT1 - UTC |
[in] | xp | "x" polar motion |
[in] | yp | "y" polar motion |
[in] | wavelength_um | Wavelength (micrometers) |
[out] | observed_ra | Apparent right ascension (returned) |
[out] | observed_dec | Apparent declination (returned) |
[out] | hour_angle | Hour angle |
Definition at line 712 of file CECoordinates.cpp.
CECoordinates CECoordinates::ConvertTo | ( | CECoordinateType | output_coord_type, |
const CEObserver & | observer, | ||
const CEDate & | date = CEDate::CurrentJD() |
||
) |
Convert these coordinates to another coordinate system NOTE: If this object is not OBSERVED coordinates, only JD is needed.
NOTE: If 'observer' is not nulltr then 'jd' will be ignored and the date will instead be derived from the observers date.
[in] | output_coord_type | Output coordinate type (see CECoordinateType) |
[in] | date | Julian date for conversion |
[in] | observer | If these coordinates are OBSERVED, then observer represents the observer for these coordinates. Otherwise they represent the observer to convert the coordinates to. |
Definition at line 1261 of file CECoordinates.cpp.
CECoordinates CECoordinates::ConvertTo | ( | CECoordinateType | output_coord_type, |
double | jd = CEDate::CurrentJD() , |
||
double | longitude = 0.0 , |
||
double | latitude = 0.0 , |
||
double | elevation_m = 0.0 , |
||
double | pressure_hPa = -1.0 , |
||
double | temperature_celsius = -1000 , |
||
double | relative_humidity = 0.0 , |
||
double | dut1 = 0.0 , |
||
double | xp = 0.0 , |
||
double | yp = 0.0 , |
||
double | wavelength_um = 0.5 |
||
) |
Convert these coordinates to another coordinate system NOTE: If neither this object or output_coord_type are OBSERVED coordinates, only JD is needed.
[in] | output_coord_type | Output coordinate type (see CECoordinateType) |
[in] | jd | Julian date for conversion |
[in] | longitude | Observer longitude (radians, east-positive) |
[in] | latitude | Observer latitude (radians) |
[in] | elevation_m | Elevation (meters above sea level) |
[in] | pressure_hPa | Atmospheric pressure (hPa) |
[in] | temperature_celsius | Atmospheric temperature (degrees celsius) |
[in] | relative_humidity | Relative humidity |
[in] | dut1 | DUT1 value representing UTC-UT1 |
[in] | xp | x-polar motion |
[in] | yp | y-polar motion |
[in] | wavelength_um | wavelength (micrometers) |
Definition at line 1300 of file CECoordinates.cpp.
CECoordinates CECoordinates::ConvertToCIRS | ( | double | jd = CEDate::CurrentJD() , |
double | longitude = 0.0 , |
||
double | latitude = 0.0 , |
||
double | elevation_m = 0.0 , |
||
double | pressure_hPa = -1.0 , |
||
double | temperature_celsius = -1000 , |
||
double | relative_humidity = 0.0 , |
||
double | dut1 = 0.0 , |
||
double | xp = 0.0 , |
||
double | yp = 0.0 , |
||
double | wavelength_um = 0.5 |
||
) |
Convert these coordinates to CIRS coordinates.
NOTE: If this object is not OBSERVED coordinates, only JD is needed.
[in] | jd | Julian date for conversion |
[in] | longitude | Observer longitude (radians, east-positive) |
[in] | latitude | Observer latitude (radians) |
[in] | elevation_m | Elevation (meters above sea level) |
[in] | pressure_hPa | Atmospheric pressure (hPa) |
[in] | temperature_celsius | Atmospheric temperature (degrees celsius) |
[in] | relative_humidity | Relative humidity |
[in] | dut1 | DUT1 value representing UTC-UT1 |
[in] | xp | x-polar motion |
[in] | yp | y-polar motion |
[in] | wavelength_um | wavelength (micrometers) |
Definition at line 1358 of file CECoordinates.cpp.
CECoordinates CECoordinates::ConvertToGalactic | ( | double | jd = CEDate::CurrentJD() , |
double | longitude = 0.0 , |
||
double | latitude = 0.0 , |
||
double | elevation_m = 0.0 , |
||
double | pressure_hPa = -1.0 , |
||
double | temperature_celsius = -1000 , |
||
double | relative_humidity = 0.0 , |
||
double | dut1 = 0.0 , |
||
double | xp = 0.0 , |
||
double | yp = 0.0 , |
||
double | wavelength_um = 0.5 |
||
) |
Convert these coordinates to GALACTIC coordinates.
NOTE: If this object is not OBSERVED coordinates, only JD is needed. NOTE: If this object is ICRS coordinates, JD is not necessary
[in] | jd | Julian date for conversion |
[in] | longitude | Observer longitude (radians, east-positive) |
[in] | latitude | Observer latitude (radians) |
[in] | elevation_m | Elevation (meters above sea level) |
[in] | pressure_hPa | Atmospheric pressure (hPa) |
[in] | temperature_celsius | Atmospheric temperature (degrees celsius) |
[in] | relative_humidity | Relative humidity |
[in] | dut1 | DUT1 value representing UTC-UT1 |
[in] | xp | x-polar motion |
[in] | yp | y-polar motion |
[in] | wavelength_um | wavelength (micrometers) |
Definition at line 1471 of file CECoordinates.cpp.
CECoordinates CECoordinates::ConvertToICRS | ( | double | jd = CEDate::CurrentJD() , |
double | longitude = 0.0 , |
||
double | latitude = 0.0 , |
||
double | elevation_m = 0.0 , |
||
double | pressure_hPa = -1.0 , |
||
double | temperature_celsius = -1000 , |
||
double | relative_humidity = 0.0 , |
||
double | dut1 = 0.0 , |
||
double | xp = 0.0 , |
||
double | yp = 0.0 , |
||
double | wavelength_um = 0.5 |
||
) |
Convert these coordinates to ICRS coordinates.
NOTE: If this object is not OBSERVED coordinates, only JD is needed. NOTE: If this object is GALACTIC coordinates, JD is not necessary
[in] | jd | Julian date for conversion |
[in] | longitude | Observer longitude (radians, east-positive) |
[in] | latitude | Observer latitude (radians) |
[in] | elevation_m | Elevation (meters above sea level) |
[in] | pressure_hPa | Atmospheric pressure (hPa) |
[in] | temperature_celsius | Atmospheric temperature (degrees celsius) |
[in] | relative_humidity | Relative humidity |
[in] | dut1 | DUT1 value representing UTC-UT1 |
[in] | xp | x-polar motion |
[in] | yp | y-polar motion |
[in] | wavelength_um | wavelength (micrometers) |
Definition at line 1418 of file CECoordinates.cpp.
CECoordinates CECoordinates::ConvertToObserved | ( | double | jd = CEDate::CurrentJD() , |
double | longitude = 0.0 , |
||
double | latitude = 0.0 , |
||
double | elevation_m = 0.0 , |
||
double | pressure_hPa = -1.0 , |
||
double | temperature_celsius = -1000 , |
||
double | relative_humidity = 0.0 , |
||
double | dut1 = 0.0 , |
||
double | xp = 0.0 , |
||
double | yp = 0.0 , |
||
double | wavelength_um = 0.5 |
||
) |
Convert these coordinates to observed coordinates.
Note: For accurate positions you MUST supply at least jd, longitude, and latitude
[in] | jd | Julian date for conversion |
[in] | longitude | Observer longitude (radians, east-positive) |
[in] | latitude | Observer latitude (radians) |
[in] | elevation_m | Elevation (meters above sea level) |
[in] | pressure_hPa | Atmospheric pressure (hPa) |
[in] | temperature_celsius | Atmospheric temperature (degrees celsius) |
[in] | relative_humidity | Relative humidity |
[in] | dut1 | DUT1 value representing UTC-UT1 |
[in] | xp | x-polar motion |
[in] | yp | y-polar motion |
[in] | wavelength_um | wavelength (micrometers) |
Definition at line 1522 of file CECoordinates.cpp.
|
private |
Copy data members from another CECoordinates object.
[in] | other | Another coordinates objec to copy |
Definition at line 1719 of file CECoordinates.cpp.
|
static |
Convert a given angle vector from {degrees, minutes, seconds} to an angle.
[in] | angle | What is the angle value |
[in] | return_type | Specifies angle type for returned value |
The supplied angle
value should be of length 3 with the indices having the following values -[0] = Degrees -[1] = Arcminutes -[2] = Arcseconds -[3] = Arcsec fraction (can be omitted if [2] is the decimal value)
Definition at line 1662 of file CECoordinates.cpp.
|
private |
Cleanup data members that need to be freed or cleared.
Definition at line 1730 of file CECoordinates.cpp.
|
static |
Galactic -> CIRS coordinate conversion.
[in] | glon | Galactic longitude |
[in] | glat | Galactic latitude |
[out] | ra | CIRS right ascension |
[out] | dec | CIRS declination |
[in] | date | Date information |
[in] | angle_type | Angle type (RADIANS or DEGREES) |
Definition at line 441 of file CECoordinates.cpp.
|
static |
Galactic -> ICRS coordinate conversion.
[in] | glon | Galactic longitude |
[in] | glat | Galactic latitude |
[out] | ra | ICRS right ascension (returned) |
[out] | dec | ICRS declinaton (returned) |
[in] | angle_type | Angle type (DEGREES or RADIANS) |
Definition at line 474 of file CECoordinates.cpp.
|
static |
Galactic -> Observed (i.e.
observer specific) coordinate conversion. For the raw version of this method see CECoordinates::Galactic2Observed().
[in] | glon | Galactic longitude |
[in] | glat | Galactic latitude |
[out] | az | Azimuth angle (returned) |
[out] | zen | Zenith angle (returned) |
[in] | observer | Observer information |
[in] | wavelength | Wavelength |
[out] | observed_glon | Observed galactic longitude (returned) |
[out] | observed_glat | Observed galactic latitude (returned) |
[out] | hour_angle | Observed hour angle of object |
Definition at line 508 of file CECoordinates.cpp.
|
static |
Raw method for converting Galactic -> Observed (observer specific) coordinates.
Note: All angles are expected to be in radians.
[in] | glon | Galactic longitude (radians) |
[in] | glat | Galactic latitude (radians) |
[out] | az | Observed azimuth angle (radians, returned) |
[out] | zen | Observed zenith angle (radians, returned) |
[in] | julian_date | Julian date for conversion |
[in] | longitude | Observer geographic longitude (radians) |
[in] | latitude | Observer geographic latitude (radians) |
[in] | elevation_m | Observer elevation (meters) |
[in] | pressure_hPa | Atmospheric pressure (HPa) |
[in] | temperature_celsius | Temperature (degrees Celsius) |
[in] | relative_humidity | Relative humidity (0.0 - 1.0) |
[in] | dut1 | UTC - UT1 |
[in] | xp | "x" polar motion |
[in] | yp | "y" polar motion |
[in] | wavelength | Wavelength |
[out] | observed_glon | Apparent galactic longitude (returned) |
[out] | observed_glat | Apparent galactic latitude (returned) |
[out] | hour_angle | Hour angle |
Definition at line 965 of file CECoordinates.cpp.
|
inline |
Return coordinate system.
Definition at line 486 of file CECoordinates.h.
|
static |
Convert a given angle into degrees, arcminutes, arcseconds.
[in] | angle | What is the angle value |
[in] | angle_type | Specifies what type of angle |
Definition at line 1594 of file CECoordinates.cpp.
|
static |
Convert a given angle into hours, minutes, seconds.
This method is almost exclusively only useful for right ascension
[in] | angle | What is the angle value |
[in] | angle_type | Specifies what type of angle |
Definition at line 1613 of file CECoordinates.cpp.
|
virtual |
Return the observed coordinates using an observer object (CEObserver)
[in] | date | Julian date of the observation |
[in] | observer | Observer information |
[in] | wavelength_um | Wavelength being observed (micrometers) |
Definition at line 1143 of file CECoordinates.cpp.
|
virtual |
Return the local sky coordinates of this object as a CECoordinates object.
[in] | julian_date | Julian date of the observation |
[in] | longitude | Observer longitude (radians) |
[in] | latitude | Observer latitude (radians) |
[in] | elevation_m | Observer elevation (meters) |
[in] | pressure_hPa | Observer atmospheric pressure (hPa) |
[in] | temerature_celsius | Temperature (degrees Celsius) |
[in] | relative_humidity | Relative humidity |
[in] | dut1 | 'UTC-UT1' |
[in] | xp | x-polar motion |
[in] | yp | y-polar motion |
[in] | wavelength_um | Wavelength being observed (micrometers) |
Definition at line 1087 of file CECoordinates.cpp.
|
static |
Convert from {hours, minutes, seconds} to an angle.
[in] | angle | What is the angle value |
[in] | return_type | Specifies angle type for returned value |
The supplied angle
value should be of length 3 with the indices having the following values -[0] = Hours -[1] = Minutes -[2] = Seconds -[3] = Fraction of a second (not absolutely necessary)
Definition at line 1636 of file CECoordinates.cpp.
|
static |
ICRS -> CIRS coordinate conversion.
Uses the SOFA function 'iauAtci13'
[in] | input_ra | Right ascension to be converted |
[in] | input_dec | Declination to be converted |
[out] | return_ra | CIRS Right ascension (returned) |
[out] | return_dec | CIRS Declination (returned) |
[in] | date | Date object |
[in] | angle_type | Angle type (either DEGREES or RADIANS) |
Definition at line 287 of file CECoordinates.cpp.
|
static |
ICRS -> Galactic coordinate conversion (uses the SOFA 'iauIcrs2g' function)
[in] | input_ra | ICRS Right ascension |
[in] | input_dec | ICRS Declination |
[out] | glon | Galactic longitude |
[out] | glat | Galactic latitude |
[in] | angle_type | Angle type |
Definition at line 336 of file CECoordinates.cpp.
|
static |
ICRS -> Observed coordinate conversion.
The integer returned is a status code with the following meanings: +1 = dubious year (too far into the past/future to be trusted) 0 = OK status -1 = unacceptable date
[in] | ra | ICRS right ascension |
[in] | dec | ICRS declination |
[out] | az | Azimuth |
[out] | zen | Zenith |
[in] | observer | Observer information |
[in] | angle_type | Angle type |
[out] | observed_ra | Apparent ICRS Right Ascension |
[out] | observed_dec | Apparent ICRS Declination |
[out] | hour_angle | Observed hour angle of object |
Definition at line 373 of file CECoordinates.cpp.
|
static |
Raw method for converting CIRS -> Observed (observer specific) coordinates.
Note: All angles are expected to be in radians.
[in] | ra | CIRS right ascension (radians) |
[in] | dec | CIRS declination (radians) |
[out] | az | Observed azimuth angle (radians, returned) |
[out] | zen | Observed zenith angle (radians, returned) |
[in] | julian_date | Julian date for conversion |
[in] | longitude | Observer geographic longitude (radians) |
[in] | latitude | Observer geographic latitude (radians) |
[in] | elevation_m | Observer elevation (meters) |
[in] | pressure_hPa | Atmospheric pressure (HPa) |
[in] | temperature_celsius | Temperature (degrees Celsius) |
[in] | relative_humidity | Relative humidity (0.0 - 1.0) |
[in] | dut1 | UTC - UT1 |
[in] | xp | "x" polar motion |
[in] | yp | "y" polar motion |
[in] | wavelength | Wavelength (micrometers) |
[out] | observed_ra | Apparent right ascension (returned) |
[out] | observed_dec | Apparent declination (returned) |
[out] | hour_angle | Hour angle |
Definition at line 838 of file CECoordinates.cpp.
|
private |
Set initial values and allocate memory for data members.
Definition at line 1739 of file CECoordinates.cpp.
|
static |
Convert Observed -> CIRS coordinates.
[in] | az | Azimuth of the coordinates |
[in] | zen | Zenith angle of the coordinates |
[out] | ra | CIRS right ascension (returned) |
[out] | dec | CIRS declination (returned) |
[in] | observer | CEObserver object describing observer |
[in] | angle_type | Angle type of input/output coordinates |
Definition at line 566 of file CECoordinates.cpp.
|
static |
Raw method for converting Observed (observer specific) -> CIRS coordinates (uses the SOFA 'iauAtoi13' function) Note: All angles are expected to be in radians.
[in] | az | Observed azimuth angle |
[in] | zen | Observed zenith angle |
[out] | ra | CIRS right ascension (returned) |
[out] | dec | CIRS declination (returned) |
[in] | julian_date | Julian date for conversion |
[in] | longitude | Observer geographic longitude |
[in] | latitude | Observer geographic latitude |
[in] | elevation_m | Observer elevation (meters) |
[in] | pressure_hPa | Atmospheric pressure (HPa) |
[in] | temperature_celsius | Temperature (degrees Celsius) |
[in] | relative_humidity | Relative humidity (0.0 - 1.0) |
[in] | dut1 | UTC - UT1 |
[in] | xp | "x" polar motion |
[in] | yp | "y" polar motion |
[in] | wavelength_um | Wavelength (micrometers) |
Definition at line 787 of file CECoordinates.cpp.
|
static |
Convert Observed -> Galactic coordinates.
[in] | az | Azimuth of the coordinates |
[in] | zen | Zenith angle of the coordinates |
[out] | glon | Galactic longitude (returned) |
[out] | glat | Galactic latitude (returned) |
[in] | observer | CEObserver object describing observer |
[in] | angle_type | Angle type of input/output coordinates |
Definition at line 653 of file CECoordinates.cpp.
|
static |
Raw method for converting Observed (observer specific) -> Galactic coordinates Note: All angles are expected to be in radians.
[in] | az | Observed azimuth angle (radians) |
[in] | zen | Observed zenith angle (radians) |
[out] | glon | Galactic longitude (radians, returned) |
[out] | glat | Galactic latitude (radians, returned) |
[in] | julian_date | Julian date for conversion |
[in] | longitude | Observer geographic longitude (radians) |
[in] | latitude | Observer geographic latitude (radians) |
[in] | elevation_m | Observer elevation (meters) |
[in] | pressure_hPa | Atmospheric pressure (HPa) |
[in] | temperature_celsius | Temperature (degrees Celsius) |
[in] | relative_humidity | Relative humidity (0.0 - 1.0) |
[in] | dut1 | UTC - UT1 |
[in] | xp | "x" polar motion |
[in] | yp | "y" polar motion |
[in] | wavelength_um | Wavelength (micrometers) |
Definition at line 1041 of file CECoordinates.cpp.
|
static |
Convert Observed -> ICRS coordinates.
[in] | az | Azimuth of the coordinates |
[in] | zen | Zenith angle of the coordinates |
[out] | ra | CIRS right ascension (returned) |
[out] | dec | CIRS declination (returned) |
[in] | observer | CEObserver object describing observer |
[in] | angle_type | Angle type of input/output coordinates |
Definition at line 608 of file CECoordinates.cpp.
|
static |
Raw method for converting Observed (observer specific) -> ICRS coordinates Note: All angles are expected to be in radians.
[in] | az | Observed azimuth angle (radians) |
[in] | zen | Observed zenith angle (radians) |
[out] | ra | ICRS right ascension (radians, returned) |
[out] | dec | ICRS declination (radians, returned) |
[in] | julian_date | Julian date for conversion |
[in] | longitude | Observer geographic longitude (radians) |
[in] | latitude | Observer geographic latitude (radians) |
[in] | elevation_m | Observer elevation (meters) |
[in] | pressure_hPa | Atmospheric pressure (HPa) |
[in] | temperature_celsius | Temperature (degrees Celsius) |
[in] | relative_humidity | Relative humidity (0.0 - 1.0) |
[in] | dut1 | UTC - UT1 |
[in] | xp | "x" polar motion |
[in] | yp | "y" polar motion |
[in] | wavelength_um | Wavelength (micrometers) |
Definition at line 911 of file CECoordinates.cpp.
CECoordinates & CECoordinates::operator= | ( | const CECoordinates & | other | ) |
Overloaded '=' (assignment) operator.
[in] | other | Coordinates object to be copied |
Definition at line 122 of file CECoordinates.cpp.
std::string CECoordinates::print | ( | void | ) | const |
Generate a message string that specifies the information about this coordinate.
Definition at line 1705 of file CECoordinates.cpp.
|
virtual |
Set the coordinates of this object.
[in] | xcoord | X-coordinate |
[in] | ycoord | Y-coordinate |
[in] | coord_type | Coordinate type (see CECoordinateType) |
[in] | angle_type | Specifies whether xcoord & ycoord are in RADIANS or DEGREES (see CEAngleType) |
Definition at line 1682 of file CECoordinates.cpp.
|
virtual |
Set the coordinates from another CECoordinates object.
[in] | coords | Another coordinates object to copy |
Definition at line 1695 of file CECoordinates.cpp.
|
inlinevirtual |
Return x coordinate at given Julian date.
[in] | jd | Julian date (used only by derived classes) |
Definition at line 411 of file CECoordinates.h.
|
inlinevirtual |
Returns x coordinate at given julian date in degrees.
[in] | jd | Julian date (used only by derived classes) |
Definition at line 450 of file CECoordinates.h.
|
inlinevirtual |
Return x coordinate at given Julian date in radians.
[in] | jd | Julian date (used only by derived classes) |
Definition at line 437 of file CECoordinates.h.
|
inlinevirtual |
Return y coordinate at given Julian date.
[in] | jd | Julian date (used only by derived classes) |
Definition at line 424 of file CECoordinates.h.
|
inlinevirtual |
Returns y coordinate at given Julian date in degrees.
[in] | jd | Julian date (used only by derived classes) |
Definition at line 475 of file CECoordinates.h.
|
inlinevirtual |
Returns y coordinate at given Julian date in radians.
[in] | jd | Julian date (used only by derived classes) |
Definition at line 463 of file CECoordinates.h.
|
friend |
Compare two coordinate objects.
Definition at line 1778 of file CECoordinates.cpp.
|
friend |
Compare two coordinate objects.
Definition at line 1751 of file CECoordinates.cpp.
|
protected |
Definition at line 387 of file CECoordinates.h.
|
mutableprotected |
Definition at line 385 of file CECoordinates.h.
|
mutableprotected |
Definition at line 386 of file CECoordinates.h.