GammaLib  2.1.0.dev
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
GSkyDir Class Reference

Sky direction class. More...

#include <GSkyDir.hpp>

Inheritance diagram for GSkyDir:
GBase

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...
 
GSkyDiroperator= (const GSkyDir &dir)
 Assignment operator. More...
 
void clear (void)
 Clear sky direction. More...
 
GSkyDirclone (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...
 

Detailed Description

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.

Constructor & Destructor Documentation

GSkyDir::GSkyDir ( void  )

Constructor.

Definition at line 60 of file GSkyDir.cpp.

References init_members().

Referenced by clone().

GSkyDir::GSkyDir ( const GVector vector)
explicit

Vector constructor.

Parameters
[in]vectorSky 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.

Parameters
[in]dirSky direction.

Definition at line 96 of file GSkyDir.cpp.

References copy_members(), and init_members().

GSkyDir::~GSkyDir ( void  )
virtual

Destructor.

Definition at line 112 of file GSkyDir.cpp.

References free_members().

Member Function Documentation

const double & GSkyDir::b ( void  ) const
inline

Return galactic latitude in radians.

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

double GSkyDir::b_deg ( void  ) const
inline

Returns galactic latitude in degrees.

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.

Parameters
[in]vector3D 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.

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

Referenced by GSkyDir(), precess(), and rotate().

std::string GSkyDir::classname ( void  ) const
inlinevirtual

Return class name.

Returns
String containing the class name ("GSkyDir").

Implements GBase.

Definition at line 148 of file GSkyDir.hpp.

GSkyDir * GSkyDir::clone ( void  ) const
virtual

Clone sky direction.

Returns
Pointer to deep copy of sky direction.

Implements GBase.

Definition at line 182 of file GSkyDir.cpp.

References GSkyDir().

void GSkyDir::copy_members ( const GSkyDir dir)
private

Copy class members.

Parameters
[in]dirSky 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.

Parameters
[in]dirSky direction to which cosine of distance is to be computed.
Returns
Cosine of angular distance.

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

const double & GSkyDir::dec ( void  ) const
inline
double GSkyDir::dist ( const GSkyDir dir) const
inline

Compute angular distance between sky directions in radians.

Parameters
[in]dirSky direction to which distance is to be computed.
Returns
Angular distance (radians).

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

void GSkyDir::equ2gal ( void  ) const
private

Convert equatorial to galactic coordinates.

Definition at line 1332 of file GSkyDir.cpp.

References b(), euler(), l(), m_b, m_dec, m_has_lb, m_l, and m_ra.

Referenced by b(), galvector(), and l().

void GSkyDir::euler ( const int &  type,
const double &  xin,
const double &  yin,
double *  xout,
double *  yout 
) const
private

General coordinate transformation routine for J2000.

Parameters
[in]typeConversion type (0=equ2gal, 1=gal2equ)
[in]xinInput longitude (RA or GLON) in radians.
[in]yinInput latitude (Dec or GLAT) in radians.
[out]xoutOutput longitude in radians.
[out]youtOutput 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.

Referenced by equ2gal(), and gal2equ().

void GSkyDir::free_members ( void  )
private

Delete class members.

Definition at line 1322 of file GSkyDir.cpp.

Referenced by clear(), operator=(), and ~GSkyDir().

void GSkyDir::gal2equ ( void  ) const
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.

Parameters
[in]vector3D 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.

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

void GSkyDir::init_members ( void  )
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=().

const double & GSkyDir::l ( void  ) const
inline

Return galactic longitude in radians.

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

double GSkyDir::l_deg ( void  ) const
inline

Return galactic longitude in degrees.

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

Parameters
[in]lGalactic longitude in radians.
[in]bGalactic 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)

Parameters
[in]lGalactic longitude in degrees.
[in]bGalactic 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.

Parameters
[in]timeTime.
[in]epochJulian 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().

GSkyDir & GSkyDir::operator= ( const GSkyDir dir)

Assignment operator.

Parameters
[in]dirSky direction.
Returns
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.

Parameters
[in]dirSky direction.
[in]coordsysCoordinate system ("CEL" or "GAL")
Returns
Position angle (rad).

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

  • \((\alpha_0,\delta_0)\) are the celestial coordinates of the instance, and
  • \((\alpha,\delta)\) are the celestial coordinates of sky direction.

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

  • \((\l_0,\l_0)\) are the Galactic coordinates of the instance, and
  • \((\l,\b)\) are the Galactic coordinates of sky direction.

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

double GSkyDir::posang_deg ( const GSkyDir dir,
const std::string &  coordsys = "CEL" 
) const
inline

Compute position angle between sky directions in degrees.

Parameters
[in]dirSky direction.
[in]coordsysCoordinate system ("CEL" or "GAL")
Returns
Position angle (deg).

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.

Parameters
[in]from_epochEpoch of the current coordinate.
[in]to_epochEpoch 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().

Referenced by moon(), and sun().

std::string GSkyDir::print ( const GChatter chatter = NORMAL) const
virtual

Print sky direction information.

Parameters
[in]chatterChattiness.
Returns
String containing sky direction information.

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

const double & GSkyDir::ra ( void  ) const
inline
void GSkyDir::radec ( const double &  ra,
const double &  dec 
)

Set equatorial sky direction (radians)

Parameters
[in]raRight Ascension in radians.
[in]decDeclination 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::rotate ( const double &  phi,
const double &  theta 
)

Rotate sky direction by zenith and azimuth angle.

Parameters
[in]phiAzimuth angle (radians).
[in]thetaZenith 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.

Parameters
[in]phiAzimuth angle (deg).
[in]thetaZenith 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.

Parameters
[in]timeTime.
[in]epochJulian 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().

Friends And Related Function Documentation

bool operator!= ( const GSkyDir a,
const GSkyDir b 
)
friend

Non equality operator.

Parameters
[in]aFirst sky direction.
[in]bSecond sky direction.

Definition at line 1495 of file GSkyDir.cpp.

bool operator== ( const GSkyDir a,
const GSkyDir b 
)
friend

Equality operator.

Parameters
[in]aFirst sky direction.
[in]bSecond 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.

Member Data Documentation

double GSkyDir::m_b
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().

double GSkyDir::m_cos_b
mutableprivate

Definition at line 135 of file GSkyDir.hpp.

Referenced by copy_members(), cos_dist(), galvector(), init_members(), and posang().

double GSkyDir::m_cos_dec
mutableprivate

Definition at line 137 of file GSkyDir.hpp.

Referenced by celvector(), copy_members(), cos_dist(), init_members(), and posang().

double GSkyDir::m_dec
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().

bool GSkyDir::m_has_lb
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().

bool GSkyDir::m_has_lb_cache
mutableprivate
bool GSkyDir::m_has_radec
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().

bool GSkyDir::m_has_radec_cache
mutableprivate
double GSkyDir::m_l
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().

double GSkyDir::m_ra
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().

double GSkyDir::m_sin_b
mutableprivate

Definition at line 134 of file GSkyDir.hpp.

Referenced by copy_members(), cos_dist(), galvector(), init_members(), and posang().

double GSkyDir::m_sin_dec
mutableprivate

Definition at line 136 of file GSkyDir.hpp.

Referenced by celvector(), copy_members(), cos_dist(), init_members(), and posang().


The documentation for this class was generated from the following files: