GammaLib
2.1.0.dev
|
Sky direction class. More...
#include <GSkyDir.hpp>
Public Member Functions | |
GSkyDir (void) | |
Constructor. More... | |
GSkyDir (const GVector &vector) | |
Vector constructor. More... | |
GSkyDir (const GSkyDir &dir) | |
Copy constructor. More... | |
virtual | ~GSkyDir (void) |
Destructor. More... | |
GSkyDir & | operator= (const GSkyDir &dir) |
Assignment operator. More... | |
void | clear (void) |
Clear sky direction. More... | |
GSkyDir * | clone (void) const |
Clone sky direction. More... | |
std::string | classname (void) const |
Return class name. More... | |
void | radec (const double &ra, const double &dec) |
Set equatorial sky direction (radians) More... | |
void | radec_deg (const double &ra, const double &dec) |
Set equatorial sky direction (degrees) More... | |
void | lb (const double &l, const double &b) |
Set galactic sky direction (radians) More... | |
void | lb_deg (const double &l, const double &b) |
Set galactic sky direction (degrees) More... | |
void | celvector (const GVector &vector) |
Set sky direction from 3D vector in celestial coordinates. More... | |
void | galvector (const GVector &vector) |
Set sky direction from 3D vector in Galactic coordinates. More... | |
void | rotate (const double &phi, const double &theta) |
Rotate sky direction by zenith and azimuth angle. More... | |
void | rotate_deg (const double &phi, const double &theta) |
Rotate sky direction by zenith and azimuth angle. More... | |
void | precess (const double &from_epoch, const double &to_epoch) |
Precess sky direction. More... | |
void | sun (const GTime &time, const double &epoch=2000.0) |
Set sky direction to direction of Sun. More... | |
void | moon (const GTime &time, const double &epoch=2000.0) |
Set sky direction to direction of Moon. More... | |
const double & | l (void) const |
Return galactic longitude in radians. More... | |
const double & | b (void) const |
Return galactic latitude in radians. More... | |
const double & | ra (void) const |
Return Right Ascension in radians. More... | |
const double & | dec (void) const |
Return Declination in radians. More... | |
double | l_deg (void) const |
Return galactic longitude in degrees. More... | |
double | b_deg (void) const |
Returns galactic latitude in degrees. More... | |
double | ra_deg (void) const |
Returns Right Ascension in degrees. More... | |
double | dec_deg (void) const |
Returns Declination in degrees. More... | |
GVector | celvector (void) const |
Return sky direction as 3D vector in celestial coordinates. More... | |
GVector | galvector (void) const |
Return sky direction as 3D vector in galactic coordinates. More... | |
double | cos_dist (const GSkyDir &dir) const |
Compute cosine of angular distance between sky directions. More... | |
double | dist (const GSkyDir &dir) const |
Compute angular distance between sky directions in radians. More... | |
double | dist_deg (const GSkyDir &dir) const |
Compute angular distance between sky directions in degrees. More... | |
double | posang (const GSkyDir &dir, const std::string &coordsys="CEL") const |
Compute position angle between sky directions in radians. More... | |
double | posang_deg (const GSkyDir &dir, const std::string &coordsys="CEL") const |
Compute position angle between sky directions in degrees. More... | |
std::string | print (const GChatter &chatter=NORMAL) const |
Print sky direction information. More... | |
Public Member Functions inherited from GBase | |
virtual | ~GBase (void) |
Destructor. More... | |
Private Member Functions | |
void | init_members (void) |
Initialise class members. More... | |
void | copy_members (const GSkyDir &dir) |
Copy class members. More... | |
void | free_members (void) |
Delete class members. More... | |
void | equ2gal (void) const |
Convert equatorial to galactic coordinates. More... | |
void | gal2equ (void) const |
Convert galactic to equatorial coordinates. More... | |
void | euler (const int &type, const double &xin, const double &yin, double *xout, double *yout) const |
General coordinate transformation routine for J2000. More... | |
Private Attributes | |
bool | m_has_lb |
Has galactic coordinates. More... | |
bool | m_has_radec |
Has equatorial coordinates. More... | |
double | m_l |
Galactic longitude in radians. More... | |
double | m_b |
Galactic latitude in radians. More... | |
double | m_ra |
Right Ascension in radians. More... | |
double | m_dec |
Declination in radians. More... | |
bool | m_has_lb_cache |
bool | m_has_radec_cache |
double | m_sin_b |
double | m_cos_b |
double | m_sin_dec |
double | m_cos_dec |
Friends | |
bool | operator== (const GSkyDir &a, const GSkyDir &b) |
Equality operator. More... | |
bool | operator!= (const GSkyDir &a, const GSkyDir &b) |
Non equality operator. More... | |
Sky direction class.
This class implements a spherical coordinate on the sky. Two coordinate systems are supported: celestial (or equatorial), defined by Right Ascension and Declination, and Galactic, defined by galactic longitude and galactic latitude. The class provides automatic conversion between both systems if required. Coordinates are stored in either of the systems (in units of radians), and conversion is performed (and stored) if requested. Coordinates can be given and returned in radians or in degrees. Note that the epoch for celestial coordinates is fixed to J2000.
Definition at line 62 of file GSkyDir.hpp.
GSkyDir::GSkyDir | ( | void | ) |
Constructor.
Definition at line 60 of file GSkyDir.cpp.
References init_members().
Referenced by clone().
|
explicit |
Vector constructor.
[in] | vector | Sky vector. |
Constructs a sky direction from a 3-element celestial vector. See the celvector() method for more information.
Definition at line 78 of file GSkyDir.cpp.
References celvector(), and init_members().
GSkyDir::GSkyDir | ( | const GSkyDir & | dir | ) |
Copy constructor.
[in] | dir | Sky direction. |
Definition at line 96 of file GSkyDir.cpp.
References copy_members(), and init_members().
|
virtual |
|
inline |
Return galactic latitude in radians.
Definition at line 187 of file GSkyDir.hpp.
References equ2gal(), m_b, m_has_lb, and m_has_radec.
Referenced by b_deg(), cos_dist(), cta_psf_radial_kerns_delta::cta_psf_radial_kerns_delta(), GHealpix::dir2pix(), equ2gal(), euler(), GHealpix::interpolator(), lb(), operator==(), posang(), and precess().
|
inline |
Returns galactic latitude in degrees.
Definition at line 202 of file GSkyDir.hpp.
References b(), and gammalib::rad2deg.
Referenced by GModelSpatialElliptical::dir(), GModelSpatialRadial::dir(), GModelSpatialPointSource::dir(), and GWcs::dir2pix().
void GSkyDir::celvector | ( | const GVector & | vector | ) |
Set sky direction from 3D vector in celestial coordinates.
[in] | vector | 3D vector. |
Convert a 3-dimensional vector in celestial coordinates into a sky direction. The transformation is given by
\[ \alpha = \arctan \left( \frac{x_1}{x_0} \right) \]
\[ \delta = \arcsin x_2 \]
Definition at line 313 of file GSkyDir.cpp.
References asin(), gammalib::atan2(), m_dec, m_has_lb, m_has_lb_cache, m_has_radec, m_has_radec_cache, and m_ra.
Referenced by com_radial_kerns_omega::eval(), com_elliptical_kerns_omega::eval(), cta_nroi_radial_kern_omega::eval(), cta_nroi_elliptical_kern_omega::eval(), cta_irf_diffuse_kern_phi::eval(), cta_nroi_diffuse_kern_phi::eval(), cta_psf_diffuse_kern_phi::eval(), GEphemerides::geo2ssb(), GCTAPointing::instdir(), GCTAPointing::skydir(), GWcs::solidangle(), GCOMBvc::tdelta(), and GCOMBvcs::tdelta().
GVector GSkyDir::celvector | ( | void | ) | const |
Return sky direction as 3D vector in celestial coordinates.
Definition at line 957 of file GSkyDir.cpp.
References cos(), gal2equ(), m_cos_dec, m_dec, m_has_lb, m_has_radec, m_has_radec_cache, m_ra, m_sin_dec, and sin().
|
inlinevirtual |
Return class name.
Implements GBase.
Definition at line 148 of file GSkyDir.hpp.
|
virtual |
Clear sky direction.
Implements GBase.
Definition at line 164 of file GSkyDir.cpp.
References free_members(), and init_members().
Referenced by GSPIEventCube::alloc_data(), GSkyMap::extract(), GCOMInstDir::init_members(), GLATInstDir::init_members(), GSPIInstDir::init_members(), GPhoton::init_members(), GLATMeanPsf::init_members(), GCTAPointing::init_members(), GModelSpatialDiffuseConst::init_members(), GPulsarEphemeris::init_members(), GCOMOad::init_members(), GWcs::init_members(), GCTAInstDir::init_members(), GModelSpatialRadial::init_members(), GModelSpatialPointSource::init_members(), GModelSpatialElliptical::init_members(), GSkyRegionRectangle::init_members(), and GSkyMap::init_members().
|
virtual |
Clone sky direction.
Implements GBase.
Definition at line 182 of file GSkyDir.cpp.
References GSkyDir().
|
private |
Copy class members.
[in] | dir | Sky direction. |
Definition at line 1294 of file GSkyDir.cpp.
References m_b, m_cos_b, m_cos_dec, m_dec, m_has_lb, m_has_lb_cache, m_has_radec, m_has_radec_cache, m_l, m_ra, m_sin_b, and m_sin_dec.
Referenced by GSkyDir(), and operator=().
double GSkyDir::cos_dist | ( | const GSkyDir & | dir | ) | const |
Compute cosine of angular distance between sky directions.
[in] | dir | Sky direction to which cosine of distance is to be computed. |
Computes the cosine of the angular distance between two sky directions in radians.
Definition at line 1027 of file GSkyDir.cpp.
References b(), cos(), dec(), l(), m_b, m_cos_b, m_cos_dec, m_dec, m_has_lb, m_has_lb_cache, m_has_radec, m_has_radec_cache, m_l, m_ra, m_sin_b, m_sin_dec, ra(), and sin().
Referenced by dist().
|
inline |
Return Declination in radians.
Definition at line 241 of file GSkyDir.hpp.
References gal2equ(), m_dec, m_has_lb, and m_has_radec.
Referenced by cos_dist(), cta_psf_radial_kerns_delta::cta_psf_radial_kerns_delta(), dec_deg(), GHealpix::dir2pix(), gal2equ(), GCOMInstDir::hash(), GLATInstDir::hash(), GSPIInstDir::hash(), GCTAInstDir::hash(), GHealpix::interpolator(), moon(), operator==(), posang(), radec(), and sun().
|
inline |
Returns Declination in degrees.
Definition at line 256 of file GSkyDir.hpp.
References dec(), and gammalib::rad2deg.
Referenced by GSkyRegionCircle::dec(), GSkyRegionRectangle::dec(), GModelSpatialElliptical::dir(), GModelSpatialRadial::dir(), GModelSpatialPointSource::dir(), GWcs::dir2pix(), GLATResponse::irf(), GCOMResponse::irf_elliptical(), GResponse::irf_elliptical(), GCOMResponse::irf_radial(), GResponse::irf_radial(), GCTAResponseIrf::nroi_elliptical(), GCTAResponseIrf::nroi_radial(), GPulsarEphemeris::print(), GSPIInstDir::print(), GPhoton::print(), GLATMeanPsf::print(), GCTAInstDir::print(), GCTAResponseCube::psf_diffuse(), GCTAPointing::update(), GCTAPointing::write(), GSkyRegionCircle::write(), GCTAEventList::write(), GSkyRegionRectangle::write(), and GCTAObservation::write_attributes().
|
inline |
Compute angular distance between sky directions in radians.
[in] | dir | Sky direction to which distance is to be computed. |
Computes the angular distance between two sky directions in radians.
Definition at line 271 of file GSkyDir.hpp.
References gammalib::acos(), and cos_dist().
Referenced by GCTAOnOffObservation::compute_arf_cut(), GCTAOnOffObservation::compute_rmf(), GModelSpatialRadialProfile::contains(), GModelSpatialRadialGauss::contains(), GModelSpatialRadialGeneralGauss::contains(), GModelSpatialRadialRing::contains(), GModelSpatialRadialDisk::contains(), GModelSpatialEllipticalDisk::contains(), GModelSpatialEllipticalGauss::contains(), GModelSpatialEllipticalGeneralGauss::contains(), GModelSpatialPointSource::contains(), GModelSpatialRadialShell::contains(), GSkyRegionRectangle::dir_to_local(), dist_deg(), GModelSpatialElliptical::eval(), GModelSpatialRadial::eval(), GCTACubeExposure::fill_cube(), GCTACubePsf::fill_cube(), GCTACubeEdisp::fill_cube(), GCTAPointing::instdir(), GCTAResponseIrf::irf(), GCTAResponseCube::irf(), GCTAResponseIrf::irf_diffuse(), GCOMResponse::irf_elliptical(), GCTAResponseCube::irf_elliptical(), GCTAResponseIrf::irf_elliptical(), GCOMResponse::irf_radial(), GCTAResponseCube::irf_radial(), GCTAResponseIrf::irf_radial(), GCTAModelAeffBackground::mc(), GCTAResponseIrf::mc(), GCTAResponseIrf::nirf(), GCTAResponseIrf::npsf(), operator!=(), operator==(), GSkyRegionMap::set_region_rectangle(), GCTABackground2D::solid_angle(), GCTABackground3D::solid_angle(), GWcs::solidangle(), and GSkyMap::solidangle().
|
inline |
Compute angular distance between sky directions in degrees.
[in] | dir | Sky direction to which distance is to be computed. |
Computes the angular distance between two sky directions in degrees.
Definition at line 286 of file GSkyDir.hpp.
References dist(), and gammalib::rad2deg.
Referenced by GCOMDri::compute_dre(), GCOMDri::cone_content(), GModelSpatialDiffuseMap::contains(), GModelSpatialDiffuseCube::contains(), GSkyRegionCircle::contains(), GModelSpatialPointSource::eval(), GCOMResponse::irf(), GCOMResponse::irf_ptsrc(), GLATResponse::irf_spatial_bin(), GLATEventCube::maxrad(), GModelSky::mc(), GModelSpatialDiffuseMap::mc_cone(), GModelSpatialDiffuseCube::mc_cone(), GModelSpatialPointSource::mc_norm(), GModelSpatialElliptical::mc_norm(), GModelSpatialRadial::mc_norm(), GSkyRegionCircle::overlaps(), GCOMEventList::read_events(), GSkyMap::region_circle(), GSkyRegionMap::set_region_circle(), and GCOMOad::theta().
|
private |
|
private |
General coordinate transformation routine for J2000.
[in] | type | Conversion type (0=equ2gal, 1=gal2equ) |
[in] | xin | Input longitude (RA or GLON) in radians. |
[in] | yin | Input latitude (Dec or GLAT) in radians. |
[out] | xout | Output longitude in radians. |
[out] | yout | Output latitude in radians. |
Definition at line 1382 of file GSkyDir.cpp.
References asin(), gammalib::atan2(), b(), cos(), gammalib::fourpi, gammalib::modulo(), sin(), and gammalib::twopi.
|
private |
Delete class members.
Definition at line 1322 of file GSkyDir.cpp.
Referenced by clear(), operator=(), and ~GSkyDir().
|
private |
Convert galactic to equatorial coordinates.
Definition at line 1354 of file GSkyDir.cpp.
References dec(), euler(), m_b, m_dec, m_has_radec, m_l, m_ra, and ra().
Referenced by celvector(), dec(), ra(), and rotate().
void GSkyDir::galvector | ( | const GVector & | vector | ) |
Set sky direction from 3D vector in Galactic coordinates.
[in] | vector | 3D vector. |
Convert a 3-dimensional vector in galactic coordinates into a sky direction. The transformation is given by
\[ l = \arctan \left( \frac{x_1}{x_0} \right) \]
\[ b = \arcsin x_2 \]
Definition at line 348 of file GSkyDir.cpp.
References asin(), gammalib::atan2(), m_b, m_has_lb, m_has_lb_cache, m_has_radec, m_has_radec_cache, and m_l.
GVector GSkyDir::galvector | ( | void | ) | const |
Return sky direction as 3D vector in galactic coordinates.
Definition at line 990 of file GSkyDir.cpp.
References cos(), equ2gal(), m_b, m_cos_b, m_has_lb, m_has_lb_cache, m_has_radec, m_l, m_sin_b, and sin().
|
private |
Initialise class members.
Definition at line 1264 of file GSkyDir.cpp.
References m_b, m_cos_b, m_cos_dec, m_dec, m_has_lb, m_has_lb_cache, m_has_radec, m_has_radec_cache, m_l, m_ra, m_sin_b, and m_sin_dec.
Referenced by clear(), GSkyDir(), and operator=().
|
inline |
Return galactic longitude in radians.
Definition at line 160 of file GSkyDir.hpp.
References equ2gal(), m_has_lb, m_has_radec, and m_l.
Referenced by cos_dist(), cta_psf_radial_kerns_delta::cta_psf_radial_kerns_delta(), GHealpix::dir2pix(), equ2gal(), GHealpix::interpolator(), l_deg(), lb(), operator==(), and posang().
|
inline |
Return galactic longitude in degrees.
Definition at line 175 of file GSkyDir.hpp.
References l(), and gammalib::rad2deg.
Referenced by GModelSpatialElliptical::dir(), GModelSpatialRadial::dir(), GModelSpatialPointSource::dir(), and GWcs::dir2pix().
void GSkyDir::lb | ( | const double & | l, |
const double & | b | ||
) |
Set galactic sky direction (radians)
[in] | l | Galactic longitude in radians. |
[in] | b | Galactic latitude in radians. |
Sets Galactic longitude and latitude in radians.
Definition at line 251 of file GSkyDir.cpp.
References b(), l(), m_b, m_has_lb, m_has_lb_cache, m_has_radec, m_has_radec_cache, and m_l.
Referenced by GHealpix::loc2dir(), GHealpix::pix2dir(), and GCOMEventList::read_events().
void GSkyDir::lb_deg | ( | const double & | l, |
const double & | b | ||
) |
Set galactic sky direction (degrees)
[in] | l | Galactic longitude in degrees. |
[in] | b | Galactic latitude in degrees. |
Sets Galactic longitude and latitude in degrees.
Definition at line 278 of file GSkyDir.cpp.
References gammalib::deg2rad, m_b, m_has_lb, m_has_lb_cache, m_has_radec, m_has_radec_cache, and m_l.
Referenced by GModelSpatialElliptical::dir(), GModelSpatialRadial::dir(), GModelSpatialPointSource::dir(), GWcs::pix2dir(), GCTAPointing::read(), GSkyRegionCircle::read(), GSkyRegionRectangle::read(), and GCOMEventList::read_events().
void GSkyDir::moon | ( | const GTime & | time, |
const double & | epoch = 2000.0 |
||
) |
Set sky direction to direction of Moon.
[in] | time | Time. |
[in] | epoch | Julian epoch. |
Sets the sky direction to the direction of the Moon at a given time
.
The equations are derived from the Chapront ELP2000/82 Lunar Theory (Chapront-Touze' and Chapront, 1983, 124, 50), as described by Jean Meeus in Chapter 47 of ``Astronomical Algorithms'' (Willmann-Bell, Richmond), 2nd edition, 1998. Meeus quotes an approximate accuracy of 10 arcsec in longitude and 4 arcsec in latitude, but he does not give the time range for this accuracy.
The method was compared to the results obtained using PyEphem for the period 2000 - 2050. Differences in Right Ascension and Declination are always smaller than 0.015 degrees, corresponding to about 1 arcmin. The differences increase over time, meaning that the largest distances are reached by the end of the period.
Definition at line 598 of file GSkyDir.cpp.
References abs(), asin(), gammalib::atan2(), cos(), dec(), gammalib::deg2rad, GTime::jd(), GTime::julian_epoch(), precess(), ra(), gammalib::rad2deg, radec_deg(), sin(), and tan().
Assignment operator.
[in] | dir | Sky direction. |
Definition at line 134 of file GSkyDir.cpp.
References copy_members(), free_members(), and init_members().
double GSkyDir::posang | ( | const GSkyDir & | dir, |
const std::string & | coordsys = "CEL" |
||
) | const |
Compute position angle between sky directions in radians.
[in] | dir | Sky direction. |
[in] | coordsys | Coordinate system ("CEL" or "GAL") |
Computes the position angle of a specified sky direction with respect to the sky direction of the instance in radians. If "CEL" is specified as coordsys
the position angle is counted counterclockwise from celestial North, if "GAL" is specified the position angle is counted counterclockwise from Galactic North.
If coordsys
is "CEL" the position angle is computed using
\[PA = \arctan \left( \frac{\sin( \alpha - \alpha_0 )} {\cos \delta_0 \tan \delta - \sin \delta_0 \cos(\alpha - \alpha_0)} \right)\]
where
If coordsys
is "GAL" the position angle is computed using
\[PA = \arctan \left( \frac{\sin( \l - \l_0 )} {\cos \b_0 \tan \b - \sin \b_0 \cos(\l - \l_0)} \right)\]
where
Definition at line 1151 of file GSkyDir.cpp.
References gammalib::atan2(), b(), cos(), dec(), G_POSANG, l(), m_cos_b, m_cos_dec, m_has_lb_cache, m_has_radec_cache, m_sin_b, m_sin_dec, ra(), sin(), and tan().
Referenced by GCTAOnOffObservation::compute_arf_cut(), GCTAOnOffObservation::compute_rmf(), GSkyRegionRectangle::dir_to_local(), GModelSpatialElliptical::eval(), GCTAPointing::instdir(), GCTAResponseCube::irf_elliptical(), GCTAResponseIrf::irf_elliptical(), GCTAResponseIrf::nroi_elliptical(), GCTAResponseIrf::nroi_radial(), posang_deg(), GSPIResponse::set_cache(), and GSkyRegionMap::set_region_rectangle().
|
inline |
Compute position angle between sky directions in degrees.
[in] | dir | Sky direction. |
[in] | coordsys | Coordinate system ("CEL" or "GAL") |
Computes the position angle of a specified sky direction with respect to the sky direction of the instance in degrees. If "CEL" is specified as coordsys
the position angle is counted counterclockwise from celestial North, if "GAL" is specified the position angle is counted counterclockwise from Galactic North.
See the posang() method for more information about the computation of the position angle.
Definition at line 309 of file GSkyDir.hpp.
References posang(), and gammalib::rad2deg.
Referenced by GCOMOad::phi().
void GSkyDir::precess | ( | const double & | from_epoch, |
const double & | to_epoch | ||
) |
Precess sky direction.
[in] | from_epoch | Epoch of the current coordinate. |
[in] | to_epoch | Epoch of the returned precessed coordinate. |
Precesses the sky direction from one epoch to another.
The method is adapted from a set of fortran subroutines based on (a) pages 30-34 of the Explanatory Supplement to the AE, (b) Lieske, et al. (1977) A&A 58, 1-16, and (c) Lieske (1979) A&A 73, 282-284.
Definition at line 451 of file GSkyDir.cpp.
References b(), celvector(), cos(), gammalib::pi, and sin().
Print sky direction information.
[in] | chatter | Chattiness. |
Implements GBase.
Definition at line 1227 of file GSkyDir.cpp.
References m_b, m_dec, m_has_lb, m_has_radec, m_l, m_ra, gammalib::rad2deg, SILENT, and gammalib::str().
Referenced by GCOMInstDir::print(), GModelSpatialDiffuseMap::print(), and GCTAEventCube::print().
|
inline |
Return Right Ascension in radians.
Definition at line 214 of file GSkyDir.hpp.
References gal2equ(), m_has_lb, m_has_radec, and m_ra.
Referenced by cos_dist(), cta_psf_radial_kerns_delta::cta_psf_radial_kerns_delta(), GHealpix::dir2pix(), gal2equ(), GCOMInstDir::hash(), GLATInstDir::hash(), GSPIInstDir::hash(), GCTAInstDir::hash(), GHealpix::interpolator(), moon(), operator==(), posang(), ra_deg(), radec(), and sun().
|
inline |
Returns Right Ascension in degrees.
Definition at line 229 of file GSkyDir.hpp.
References ra(), and gammalib::rad2deg.
Referenced by GModelSpatialElliptical::dir(), GModelSpatialRadial::dir(), GModelSpatialPointSource::dir(), GWcs::dir2pix(), GLATResponse::irf(), GCOMResponse::irf_elliptical(), GResponse::irf_elliptical(), GCOMResponse::irf_radial(), GResponse::irf_radial(), GCTAResponseIrf::nroi_elliptical(), GCTAResponseIrf::nroi_radial(), GPulsarEphemeris::print(), GSPIInstDir::print(), GPhoton::print(), GLATMeanPsf::print(), GCTAInstDir::print(), GCTAResponseCube::psf_diffuse(), GSkyRegionCircle::ra(), GSkyRegionRectangle::ra(), GCTAPointing::update(), GCTAPointing::write(), GSkyRegionCircle::write(), GCTAEventList::write(), GSkyRegionRectangle::write(), and GCTAObservation::write_attributes().
void GSkyDir::radec | ( | const double & | ra, |
const double & | dec | ||
) |
Set equatorial sky direction (radians)
[in] | ra | Right Ascension in radians. |
[in] | dec | Declination in radians. |
Sets Right Ascension and Declination in radians.
Definition at line 197 of file GSkyDir.cpp.
References dec(), m_dec, m_has_lb, m_has_lb_cache, m_has_radec, m_has_radec_cache, m_ra, and ra().
Referenced by GHealpix::loc2dir(), GHealpix::pix2dir(), GCOMOads::read(), GCTABackground3D::solid_angle(), and GCTABackground2D::solid_angle().
void GSkyDir::radec_deg | ( | const double & | ra, |
const double & | dec | ||
) |
Set equatorial sky direction (degrees)
[in] | ra | Right Ascension in degrees. |
[in] | dec | Declination in degrees. |
Sets Right Ascension and Declination in degrees.
Definition at line 224 of file GSkyDir.cpp.
References gammalib::deg2rad, m_dec, m_has_lb, m_has_lb_cache, m_has_radec, m_has_radec_cache, and m_ra.
Referenced by GSkyRegionCircle::centre(), GSkyRegionRectangle::centre(), GCOMDri::compute_dre(), GCOMDri::compute_drg(), GCOMDris::compute_drws_energy(), GCOMDris::compute_drws_phibar(), GModelSpatialElliptical::dir(), GModelSpatialRadial::dir(), GModelSpatialPointSource::dir(), GPulsar::load_fermi(), GPulsar::load_integral(), GPulsar::load_parfile(), GPulsar::load_psrtime(), moon(), GWcs::pix2dir(), GCTAPointing::read(), GCTAEventList::read(), GSkyRegionCircle::read(), GSkyRegionRectangle::read(), GCTAObservation::read_attributes(), GCTAEventCube::read_cntmap(), gammalib::read_ds_roi(), GCTAEventList::read_events(), GSPIEventCube::read_pnt(), sun(), GCOMDris::vetorate_generate(), and GCOMDris::vetorate_setup().
void GSkyDir::rotate | ( | const double & | phi, |
const double & | theta | ||
) |
Rotate sky direction by zenith and azimuth angle.
[in] | phi | Azimuth angle (radians). |
[in] | theta | Zenith angle (radians). |
Rotate sky direction by a zenith and azimuth angle given in the system of the sky direction and aligned in celestial coordinates. The azimuth angle is counted counter clockwise from celestial north (this is identical to the astronomical definition of a position angle).
Definition at line 378 of file GSkyDir.cpp.
References celvector(), cos(), GMatrix::eulery(), GMatrix::eulerz(), gal2equ(), m_dec, m_has_lb, m_has_radec, m_ra, gammalib::rad2deg, and sin().
Referenced by GModelSpatial::circle_int_kern_omega::eval(), GCOMResponse::irf_diffuse(), GSkyRegionRectangle::local_to_dir(), and rotate_deg().
void GSkyDir::rotate_deg | ( | const double & | phi, |
const double & | theta | ||
) |
Rotate sky direction by zenith and azimuth angle.
[in] | phi | Azimuth angle (deg). |
[in] | theta | Zenith angle (deg). |
Rotate sky direction by a zenith and azimuth angle given in the system of the sky direction and aligned in celestial coordinates. The azimuth angle is counted counter clockwise from celestial north (this is identical to the astronomical definition of a position angle).
Definition at line 424 of file GSkyDir.cpp.
References gammalib::deg2rad, and rotate().
Referenced by GSkyRegionMap::contains(), GModelSpatialDiffuseConst::mc(), GModelSpatialRadialProfile::mc(), GCTAModelAeffBackground::mc(), GModelSpatialDiffuseMap::mc(), GModelSpatialRadialGauss::mc(), GModelSpatialRadialDisk::mc(), GModelSpatialRadialGeneralGauss::mc(), GModelSpatialRadialRing::mc(), GModelSpatialEllipticalDisk::mc(), GModelSpatialEllipticalGauss::mc(), GModelSpatialEllipticalGeneralGauss::mc(), GModelSpatialDiffuseCube::mc(), GModelSpatialRadialShell::mc(), and GCTAResponseIrf::mc().
void GSkyDir::sun | ( | const GTime & | time, |
const double & | epoch = 2000.0 |
||
) |
Set sky direction to direction of Sun.
[in] | time | Time. |
[in] | epoch | Julian epoch. |
Sets the sky direction to the direction of the Sun at a given time
.
The equations were taken from the Astronomical Almanac. They calculate the apparent coordinates of the Sun to a precision of about 0.01 degrees (36 arcsec), for dates between 1950 and 2050.
See https://en.wikipedia.org/wiki/Position_of_the_Sun
Definition at line 528 of file GSkyDir.cpp.
References asin(), gammalib::atan2(), cos(), dec(), gammalib::deg2rad, GTime::jd(), GTime::julian_epoch(), precess(), ra(), gammalib::rad2deg, radec_deg(), and sin().
Non equality operator.
[in] | a | First sky direction. |
[in] | b | Second sky direction. |
Definition at line 1495 of file GSkyDir.cpp.
Equality operator.
[in] | a | First sky direction. |
[in] | b | Second sky direction. |
Compare two sky directions. If the coordinate is at the pole, the Right Ascension or Longitude value is irrelevant.
Comparisons are done dependent on the available coordinate system. This speeds up things and avoids unnecessary coordinate transformations.
Definition at line 1431 of file GSkyDir.cpp.
|
private |
Galactic latitude in radians.
Definition at line 126 of file GSkyDir.hpp.
Referenced by b(), copy_members(), cos_dist(), equ2gal(), gal2equ(), galvector(), init_members(), lb(), lb_deg(), operator==(), and print().
|
mutableprivate |
Definition at line 135 of file GSkyDir.hpp.
Referenced by copy_members(), cos_dist(), galvector(), init_members(), and posang().
|
mutableprivate |
Definition at line 137 of file GSkyDir.hpp.
Referenced by celvector(), copy_members(), cos_dist(), init_members(), and posang().
|
private |
Declination in radians.
Definition at line 128 of file GSkyDir.hpp.
Referenced by celvector(), copy_members(), cos_dist(), dec(), equ2gal(), gal2equ(), init_members(), operator==(), print(), radec(), radec_deg(), and rotate().
|
mutableprivate |
Has galactic coordinates.
Definition at line 123 of file GSkyDir.hpp.
Referenced by b(), celvector(), copy_members(), cos_dist(), dec(), equ2gal(), galvector(), init_members(), l(), lb(), lb_deg(), operator==(), print(), ra(), radec(), radec_deg(), and rotate().
|
mutableprivate |
Definition at line 132 of file GSkyDir.hpp.
Referenced by celvector(), copy_members(), cos_dist(), galvector(), init_members(), lb(), lb_deg(), posang(), radec(), and radec_deg().
|
mutableprivate |
Has equatorial coordinates.
Definition at line 124 of file GSkyDir.hpp.
Referenced by b(), celvector(), copy_members(), cos_dist(), dec(), gal2equ(), galvector(), init_members(), l(), lb(), lb_deg(), operator==(), print(), ra(), radec(), radec_deg(), and rotate().
|
mutableprivate |
Definition at line 133 of file GSkyDir.hpp.
Referenced by celvector(), copy_members(), cos_dist(), galvector(), init_members(), lb(), lb_deg(), posang(), radec(), and radec_deg().
|
private |
Galactic longitude in radians.
Definition at line 125 of file GSkyDir.hpp.
Referenced by copy_members(), cos_dist(), equ2gal(), gal2equ(), galvector(), init_members(), l(), lb(), lb_deg(), operator==(), and print().
|
private |
Right Ascension in radians.
Definition at line 127 of file GSkyDir.hpp.
Referenced by celvector(), copy_members(), cos_dist(), equ2gal(), gal2equ(), init_members(), operator==(), print(), ra(), radec(), radec_deg(), and rotate().
|
mutableprivate |
Definition at line 134 of file GSkyDir.hpp.
Referenced by copy_members(), cos_dist(), galvector(), init_members(), and posang().
|
mutableprivate |
Definition at line 136 of file GSkyDir.hpp.
Referenced by celvector(), copy_members(), cos_dist(), init_members(), and posang().