libTheSky
Routines to compute sky positions of Sun, Moon, planets and more
All Namespaces Files Functions Variables Pages
thesky_coordinates Module Reference

Procedures for coordinates. More...

Functions/Subroutines

subroutine hc_spher_2_gc_rect (l, b, r, l0, b0, r0, x, y, z)
 Compute the geocentric rectangular coordinates of a planet, from its and the Earth's heliocentric spherical position.
 
subroutine ecl_spher_2_eq_rect (l, b, r, eps, x, y, z)
 Convert spherical, ecliptical coordinates to rectangular, equatorial coordinates of an object - both geocentric.
 
subroutine calcsunxyz (t1, l0, b0, r0, x, y, z)
 Compute the geocentric equatorial rectangular coordinates of the Sun, from Earth's heliocentric, spherical position.
 
subroutine precess_xyz (jd1, jd2, x, y, z)
 Compute the precession of the equinoxes in rectangular coordinates, from jd1 to jd2.
 
subroutine precess_eq (jd1, jd2, a1, d1)
 Compute the precession of the equinoxes in equatorial coordinates, from jd1 to jd2.
 
subroutine precess_eq_yr (yr1, yr2, a1, d1)
 Compute the precession of the equinoxes in equatorial coordinates, from yr1 to yr2.
 
subroutine precess_ecl (jd1, jd2, l, b)
 Compute the precession of the equinoxes in geocentric ecliptical coordinates.
 
subroutine precess_orb (jd1, jd2, i, o1, o2)
 Compute the precession of the equinoxes in orbital elements.
 
subroutine rect_2_spher (x, y, z, l, b, r)
 Convert rectangular coordinates x,y,z to spherical coordinates l,b,r.
 
subroutine aberration_ecl (t, l0, l, b)
 Correct ecliptic longitude and latitiude for annual aberration.
 
subroutine aberration_eq (jd, ra, dec, dra, ddec, eps0)
 Correct equatorial coordinates for annual aberration - moderate accuracy, use for stars.
 
subroutine fk5 (t, l, b)
 Convert coordinates to the FK5 system.
 
subroutine ecl_2_eq (l, b, eps, ra, dec)
 Convert (geocentric) spherical ecliptical coordinates l,b (and eps) to spherical equatorial coordinates RA, Dec.
 
subroutine eq_2_ecl (ra, dec, eps, l, b)
 Convert (geocentric) spherical equatorial coordinates RA, Dec (and eps) to spherical ecliptical coordinates l,b.
 
subroutine eq2horiz (ra, dec, agst, hh, az, alt, lat, lon)
 Convert spherical equatorial coordinates (RA, dec, agst) to spherical horizontal coordinates (hh, az, alt)
 
subroutine horiz2eq (az, alt, agst, hh, ra, dec, lat, lon)
 Convert spherical horizontal coordinates (az, alt, agst) to spherical equatorial coordinates (hh, RA, dec)
 
subroutine eq2gal (ra, dec, l, b)
 Convert spherical equatorial coordinates (RA, dec) to spherical galactic coordinates (l,b), for J2000.0!!!
 
subroutine gal2eq (l, b, ra, dec)
 Convert spherical galactic coordinates (l,b) to spherical equatorial coordinates (RA, dec), for J2000.0!!!
 
subroutine geoc2topoc_ecl (gcl, gcb, gcr, gcs, eps, lst, tcl, tcb, tcs, lat, hgt)
 Convert spherical ecliptical coordinates from the geocentric to the topocentric system.
 
subroutine geoc2topoc_eq (gcra, gcd, gcr, gch, tcra, tcd, lat, hgt)
 Convert geocentric equatorial coordinates to topocentric.
 
real(double) function refract (alt, press, temp)
 Compute the atmospheric refraction for a given true altitude. You should add the result to the uncorrected altitude in order to obtain the observed altitude. Return 0 if alt + refract < -0.3°.
 
real function refract_sp (alt, press, temp)
 Compute the atmospheric refraction for a given true altitude, using single-precision values. This is a wrapper for refract().
 
real(double) function atmospheric_refraction (alt0, h0, lat0, t0, p0, rh, lam, dtdh, eps)
 Compute the atmospheric refraction of light for a given true altitude. Return 0 for alt<-0.9°. This is a wrapper for aref(), which does the opposite (compute refraction for an apparent zenith angle). This is an expensive way to go about(!)
 
real(double) function atmospheric_refraction_apparent (alt0, h0, ph, t0, p0, rh, lam, dtdh, eps)
 Compute the atmospheric refraction of light for a given apparent (observed) altitude of a celestial object and a given observer. Note that you will usually want to use the true altitude as input instead. The method is based on N.A.O Technical Notes 59 and 63 and a paper by Standish and Auer 'Astronomical Refraction: Computational Method for all Zenith Angles'. Return 0 if alt<-1.102°.
 
real(double) function aref (z0, h0, ph, t0, p0, rh, lam, dtdh, eps)
 Compute the atmospheric refraction of light for a given apparent (observed) zenith angle. The method is based on N.A.O Technical Notes 59 and 63 and a paper by Standish and Auer 'Astronomical Refraction: Computational Method for all Zenith Angles'. Return 0 if z0>91.102°.
 
real(double) function refi (r, n, dndr)
 The refraction integrand for atmospheric_refraction_apparent()
 
subroutine troposphere_model (r0, t0, a, r, t, n, dndr)
 Troposphere model for atmospheric_refraction_apparent()
 
subroutine stratosphere_model (rt, tt, nt, a, r, n, dndr)
 Stratosphere model for atmospheric_refraction_apparent()
 

Detailed Description

Procedures for coordinates.

Function/Subroutine Documentation

◆ aberration_ecl()

subroutine thesky_coordinates::aberration_ecl ( real(double), intent(in) t,
real(double), intent(in) l0,
real(double), intent(inout) l,
real(double), intent(inout) b )

Correct ecliptic longitude and latitiude for annual aberration.

Parameters
tDynamical time in Julian Millennia since 2000.0
l0Earth longitude (rad)
lLongitude of the object (rad, I/O) (output)
bLatitude of the object (rad, I/O) (output)
See also
Meeus, Astronomical Algorithms, 1998, Ch.23

Definition at line 426 of file coordinates.f90.

Referenced by thesky_planets::planet_position().

◆ aberration_eq()

subroutine thesky_coordinates::aberration_eq ( real(double), intent(in) jd,
real(double), intent(in) ra,
real(double), intent(in) dec,
real(double), intent(out) dra,
real(double), intent(out) ddec,
real(double), intent(in), optional eps0 )

Correct equatorial coordinates for annual aberration - moderate accuracy, use for stars.

Parameters
jdJulian day of epoch
raRight ascension of the object (rad)
decDeclination of the object (rad)
draCorrection in right ascension due to aberration (rad) (output)
ddecCorrection in declination due to aberration(rad) (output)
eps0The mean obliquity of the ecliptic
See also
Meeus, Astronomical Algorithms, 1998, Ch.23

Definition at line 470 of file coordinates.f90.

References thesky_nutation::nutation().

Referenced by thesky_stars::calcstars().

Here is the call graph for this function:

◆ aref()

real(double) function thesky_coordinates::aref ( real(double), intent(in) z0,
real(double), intent(in) h0,
real(double), intent(in) ph,
real(double), intent(in) t0,
real(double), intent(in) p0,
real(double), intent(in) rh,
real(double), intent(in) lam,
real(double), intent(in) dtdh,
real(double), intent(in) eps )

Compute the atmospheric refraction of light for a given apparent (observed) zenith angle. The method is based on N.A.O Technical Notes 59 and 63 and a paper by Standish and Auer 'Astronomical Refraction: Computational Method for all Zenith Angles'. Return 0 if z0>91.102°.

Parameters
z0The observed zenith distance of the object in degrees(!)
h0The height of the observer above sea level in metres
phThe latitude of the observer in degrees(!)
t0The temperature at the observer in Kelvin
p0The pressure at the observer in millibars
rhThe relative humidity at the observer as a fraction (0-1)
lamThe wavelength of the light at the observer in micrometres
dTdhThe temperature lapse rate dT/dh in Kelvin/metre in the troposphere; the absolute value is used
epsThe desired precision in arcseconds(!)
Return values
arefThe refraction at the observer in degrees(!)
See also
Hohenkerk & Sinclair, HMNAO technical note 63 (1985)

Definition at line 1277 of file coordinates.f90.

References aref(), and atmospheric_refraction_apparent().

Referenced by aref(), and atmospheric_refraction().

Here is the call graph for this function:

◆ atmospheric_refraction()

real(double) function thesky_coordinates::atmospheric_refraction ( real(double), intent(in) alt0,
real(double), intent(in) h0,
real(double), intent(in) lat0,
real(double), intent(in) t0,
real(double), intent(in) p0,
real(double), intent(in) rh,
real(double), intent(in) lam,
real(double), intent(in) dtdh,
real(double), intent(in) eps )

Compute the atmospheric refraction of light for a given true altitude. Return 0 for alt<-0.9°. This is a wrapper for aref(), which does the opposite (compute refraction for an apparent zenith angle). This is an expensive way to go about(!)

Parameters
alt0The true (theoretical, computed) altitude of the object in radians
h0The height of the observer above sea level in metres
lat0The latitude of the observer in radians
t0The temperature at the observer in degrees Celsius (e.g. 10)
p0The pressure at the observer in millibars (e.g. 1010)
rhThe relative humidity at the observer in percent (e.g. 50)
lamThe wavelength of the light at the observer in nanometres (e.g. 550)
dTdhThe temperature lapse rate dT/dh in Kelvin/metre in the troposphere (only the absolute value is used; e.g. 0.0065)
epsThe desired precision in arcseconds (e.g. 1.d-3)
Return values
atmospheric_refractionThe refraction at the observer in radians
Todo
Adapt aref() to compute the integral in the other direction for a direct method(?)

Definition at line 1003 of file coordinates.f90.

References aref(), and atmospheric_refraction().

Referenced by atmospheric_refraction().

Here is the call graph for this function:

◆ atmospheric_refraction_apparent()

real(double) function thesky_coordinates::atmospheric_refraction_apparent ( real(double), intent(in) alt0,
real(double), intent(in) h0,
real(double), intent(in) ph,
real(double), intent(in), optional t0,
real(double), intent(in), optional p0,
real(double), intent(in), optional rh,
real(double), intent(in), optional lam,
real(double), intent(in), optional dtdh,
real(double), intent(in), optional eps )

Compute the atmospheric refraction of light for a given apparent (observed) altitude of a celestial object and a given observer. Note that you will usually want to use the true altitude as input instead. The method is based on N.A.O Technical Notes 59 and 63 and a paper by Standish and Auer 'Astronomical Refraction: Computational Method for all Zenith Angles'. Return 0 if alt<-1.102°.

Parameters
alt0The observed (apparent) altitude of the object (rad).
h0The height of the observer above sea level (metres).
phThe latitude of the observer (rad).
t0The temperature at the observer's site; optional, default 10 (°C).
p0The air pressure at the observer's site; optional, default: 1010 (mbar).
rhThe relative humidity at the observer's site; optional, default: 0.5 (fraction: 0-1).
lamThe wavelength of the light; optional, default: 550 (nm).
dTdhThe temperature gradient dT/dh in the troposphere; the absolute value is used; optional, default: 6.5e-3 (Kelvin/metre).
epsThe desired precision; optional, default: 1e-4; don't make this smaller than ~1e-8 (rad).
Return values
atmospheric_refraction_apparentThe refraction at the observer (rad).
See also
Hohenkerk & Sinclair, HMNAO technical note 63 (1985).

Definition at line 1064 of file coordinates.f90.

References atmospheric_refraction_apparent(), refi(), stratosphere_model(), and troposphere_model().

Referenced by aref(), and atmospheric_refraction_apparent().

Here is the call graph for this function:

◆ calcsunxyz()

subroutine thesky_coordinates::calcsunxyz ( real(double), intent(in) t1,
real(double), intent(in) l0,
real(double), intent(in) b0,
real(double), intent(in) r0,
real(double), intent(out) x,
real(double), intent(out) y,
real(double), intent(out) z )

Compute the geocentric equatorial rectangular coordinates of the Sun, from Earth's heliocentric, spherical position.

Parameters
t1Dynamical time in Julian Millennia after 2000.0
l0Heliocentric longitude of Earth
b0Heliocentric latitude of Earth
r0Heliocentric distance of Earth
xGeocentric rectangular X of planet (output)
yGeocentric rectangular Y of planet (output)
zGeocentric rectangular Z of planet (output)
Todo
To be disposed, has been replcaed by ecl_spher_2_eq_rect above

Definition at line 103 of file coordinates.f90.

References fk5(), and thesky_nutation::nutation().

Referenced by thesky_asteroids::asteroid_eq().

Here is the call graph for this function:

◆ ecl_2_eq()

subroutine thesky_coordinates::ecl_2_eq ( real(double), intent(in) l,
real(double), intent(in) b,
real(double), intent(in) eps,
real(double), intent(out) ra,
real(double), intent(out) dec )

Convert (geocentric) spherical ecliptical coordinates l,b (and eps) to spherical equatorial coordinates RA, Dec.

Parameters
lLongitude (rad)
bLatitude (rad)
epsObliquity of the ecliptic
raRight ascension (rad) (output)
decDeclination (rad) (output)
See also
Expl. Supl. t.t. Astronimical Almanac 3rd Ed, Eq.14.43

Definition at line 571 of file coordinates.f90.

Referenced by thesky_visibility::comet_invisible(), thesky_planets::jupiterphys(), thesky_moon::moonpos_la(), thesky_planets::planet_position(), and thesky_planets::saturnphys().

◆ ecl_spher_2_eq_rect()

subroutine thesky_coordinates::ecl_spher_2_eq_rect ( real(double), intent(in) l,
real(double), intent(in) b,
real(double), intent(in) r,
real(double), intent(in) eps,
real(double), intent(out) x,
real(double), intent(out) y,
real(double), intent(out) z )

Convert spherical, ecliptical coordinates to rectangular, equatorial coordinates of an object - both geocentric.

Parameters
lHeliocentric longitude of planet
bHeliocentric latitude of planet
rHeliocentric distance of planet
epsObliquity of the ecliptic
xGeocentric rectangular X of planet (output)
yGeocentric rectangular Y of planet (output)
zGeocentric rectangular Z of planet (output)

Definition at line 73 of file coordinates.f90.

Referenced by thesky_comets::cometgc().

◆ eq2gal()

subroutine thesky_coordinates::eq2gal ( real(double), intent(in) ra,
real(double), intent(in) dec,
real(double), intent(out) l,
real(double), intent(out) b )

Convert spherical equatorial coordinates (RA, dec) to spherical galactic coordinates (l,b), for J2000.0!!!

Parameters
raRight ascension
decDeclination
lLongitude (output)
bLatitude (output)
See also
Meeus, Astronomical Algorithms, 1998, Ch.13

Definition at line 724 of file coordinates.f90.

◆ eq2horiz()

subroutine thesky_coordinates::eq2horiz ( real(double), intent(in) ra,
real(double), intent(in) dec,
real(double), intent(in) agst,
real(double), intent(out) hh,
real(double), intent(out) az,
real(double), intent(out) alt,
real(double), intent(in), optional lat,
real(double), intent(in), optional lon )

Convert spherical equatorial coordinates (RA, dec, agst) to spherical horizontal coordinates (hh, az, alt)

Parameters
raRight ascension
decDeclination
agstGreenwich sidereal time
hhLocal hour angle (output)
azAzimuth (output)
altAltitude (output)
latGeographical latitude (rad), optional
lonGeographical longitude (rad), optional
See also
Expl. Supl. t.t. Astronimical Almanac 3rd Ed, Eq.14.44, 14.45

Definition at line 629 of file coordinates.f90.

References thesky_local::lat0, and thesky_local::lon0.

Referenced by thesky_moon::moonpos_la(), thesky_planets::planet_position(), thesky_planets::planet_position_la(), and thesky_sun::sunpos_la().

◆ eq_2_ecl()

subroutine thesky_coordinates::eq_2_ecl ( real(double), intent(in) ra,
real(double), intent(in) dec,
real(double), intent(in) eps,
real(double), intent(out) l,
real(double), intent(out) b )

Convert (geocentric) spherical equatorial coordinates RA, Dec (and eps) to spherical ecliptical coordinates l,b.

Parameters
raRight ascension (rad)
decDeclination (rad)
epsObliquity of the ecliptic
lLongitude (rad) (output)
bLatitude (rad) (output)
See also
Expl. Supl. t.t. Astronimical Almanac 3rd Ed, Eq.14.42

Definition at line 598 of file coordinates.f90.

Referenced by thesky_stars::calcstars(), thesky_comets::cometgc(), and thesky_planets::jupiterphys().

◆ fk5()

subroutine thesky_coordinates::fk5 ( real(double), intent(in) t,
real(double), intent(inout) l,
real(double), intent(inout) b )

Convert coordinates to the FK5 system.

Parameters
tDynamical time in Julian Millennia since 2000.0
lLongitude (rad, I/O) (output)
bLatitude (rad, I/O) (output)

Definition at line 538 of file coordinates.f90.

Referenced by calcsunxyz(), thesky_comets::cometgc(), and thesky_planets::planet_position().

◆ gal2eq()

subroutine thesky_coordinates::gal2eq ( real(double), intent(in) l,
real(double), intent(in) b,
real(double), intent(out) ra,
real(double), intent(out) dec )

Convert spherical galactic coordinates (l,b) to spherical equatorial coordinates (RA, dec), for J2000.0!!!

Parameters
lLongitude
bLatitude
raRight ascension (output)
decDeclination (output)

Definition at line 753 of file coordinates.f90.

◆ geoc2topoc_ecl()

subroutine thesky_coordinates::geoc2topoc_ecl ( real(double), intent(in) gcl,
real(double), intent(in) gcb,
real(double), intent(in) gcr,
real(double), intent(in) gcs,
real(double), intent(in) eps,
real(double), intent(in) lst,
real(double), intent(out) tcl,
real(double), intent(out) tcb,
real(double), intent(out) tcs,
real(double), intent(in), optional lat,
real(double), intent(in), optional hgt )

Convert spherical ecliptical coordinates from the geocentric to the topocentric system.

Parameters
gclGeocentric ecliptic longitude of the object (rad)
gcbGeocentric ecliptic latitude of the object (rad)
gcrGeocentric distance of the object
gcsGeocentric semi-diameter of the object (rad)
epsObliquity of the ecliptic (rad)
lstLocal sidereal time (rad)
tclTopocentric ecliptic longitude of the object (rad) (output)
tcbTopocentric ecliptic latitude of the object (rad) (output)
tcsTopocentric semi-diameter of the object (rad) (output)
latLatitude of the observer (rad, optional)
hgtAltitude/elevation of the observer above sea level (metres, optional)
Note
lat/lat0 and hgt/height can be provided through the module TheSky_local (rad, and m), or through the optional arguments. Note that using the latter will update the former!
See also
Meeus, Astronomical Algorithms, 1998, Ch. 11 and 40

Definition at line 797 of file coordinates.f90.

References thesky_local::height, and thesky_local::lat0.

Referenced by thesky_planets::planet_position(), and thesky_planets::planet_position_la().

◆ geoc2topoc_eq()

subroutine thesky_coordinates::geoc2topoc_eq ( real(double), intent(in) gcra,
real(double), intent(in) gcd,
real(double), intent(in) gcr,
real(double), intent(in) gch,
real(double), intent(out) tcra,
real(double), intent(out) tcd,
real(double), intent(in), optional lat,
real(double), intent(in), optional hgt )

Convert geocentric equatorial coordinates to topocentric.

Parameters
gcraGeocentric right ascension
gcdGeocentric declination
gcrGeocentric distance
gchGeocentric hour angle
tcraTopocentric right ascension (output)
tcdTopocentric declination (output)
latLatitude of the observer (rad, optional)
hgtAltitude/elevation of the observer above sea level (metres, optional)
Note
lat/lat0 and hgt/height can be provided through the module TheSky_local (rad, and m), or through the optional arguments. Note that using the latter will update the former!
See also
Meeus, Astronomical Algorithms, 1998, Ch. 11 and 40

Definition at line 858 of file coordinates.f90.

References thesky_local::height, and thesky_local::lat0.

Referenced by thesky_planets::planet_position_la().

◆ hc_spher_2_gc_rect()

subroutine thesky_coordinates::hc_spher_2_gc_rect ( real(double), intent(in) l,
real(double), intent(in) b,
real(double), intent(in) r,
real(double), intent(in) l0,
real(double), intent(in) b0,
real(double), intent(in) r0,
real(double), intent(out) x,
real(double), intent(out) y,
real(double), intent(out) z )

Compute the geocentric rectangular coordinates of a planet, from its and the Earth's heliocentric spherical position.

Parameters
lHeliocentric longitude of planet
bHeliocentric latitude of planet
rHeliocentric distance of planet
l0Heliocentric longitude of Earth
b0Heliocentric latitude of Earth
r0Heliocentric distance of Earth
xGeocentric rectangular X of planet (output)
yGeocentric rectangular Y of planet (output)
zGeocentric rectangular Z of planet (output)

Definition at line 46 of file coordinates.f90.

Referenced by thesky_planets::planet_position(), and thesky_planets::saturnphys().

◆ horiz2eq()

subroutine thesky_coordinates::horiz2eq ( real(double), intent(in) az,
real(double), intent(in) alt,
real(double), intent(in) agst,
real(double), intent(out) hh,
real(double), intent(out) ra,
real(double), intent(out) dec,
real(double), intent(in), optional lat,
real(double), intent(in), optional lon )

Convert spherical horizontal coordinates (az, alt, agst) to spherical equatorial coordinates (hh, RA, dec)

Parameters
azAzimuth
altAltitude
agstGreenwich sidereal time
hhLocal hour angle (output)
raRight ascension (output)
decDeclination (output)
latGeographical latitude (rad), optional
lonGeographical longitude (rad), optional
See also
Expl. Supl. t.t. Astronimical Almanac 3rd Ed, Eq.14.46, 14.47

Definition at line 680 of file coordinates.f90.

References thesky_local::lat0, and thesky_local::lon0.

Referenced by thesky_visibility::limmag_zenith_jd().

◆ precess_ecl()

subroutine thesky_coordinates::precess_ecl ( real(double), intent(in) jd1,
real(double), intent(in) jd2,
real(double), intent(inout) l,
real(double), intent(inout) b )

Compute the precession of the equinoxes in geocentric ecliptical coordinates.

Parameters
jd1Original Julian day
jd2Target Julian day
lEcliptic longitude (I/O, rad) (output)
bEcliptic latitude (I/O, rad) (output)
See also
Meeus, Astronomical Algorithms, 1998, Ch.21, p.136

Definition at line 295 of file coordinates.f90.

Referenced by thesky_moon::elp_mpp02_lbr(), thesky_moonroutines::moonphys(), and thesky_planets::planet_position().

◆ precess_eq()

subroutine thesky_coordinates::precess_eq ( real(double), intent(in) jd1,
real(double), intent(in) jd2,
real(double), intent(inout) a1,
real(double), intent(inout) d1 )

Compute the precession of the equinoxes in equatorial coordinates, from jd1 to jd2.

Parameters
jd1Original Julian day
jd2Target Julian day
a1Right ascension (IO, rad) (output)
d1Declination (IO, rad) (output)
See also
Meeus, Astronomical Algorithms, 1998, Ch.21, p.134

Definition at line 210 of file coordinates.f90.

Referenced by thesky_stars::calcstars(), thesky_functions::const_id(), and thesky_moonroutines::moonphys().

◆ precess_eq_yr()

subroutine thesky_coordinates::precess_eq_yr ( real(double), intent(in) yr1,
real(double), intent(in) yr2,
real(double), intent(inout) a1,
real(double), intent(inout) d1 )

Compute the precession of the equinoxes in equatorial coordinates, from yr1 to yr2.

Parameters
yr1Original year
yr2Target year
a1Right ascension (IO, rad) (output)
d1Declination (IO, rad) (output)
See also
Meeus, Astronomical Algorithms, 1998, Ch.21, p.134

Definition at line 253 of file coordinates.f90.

◆ precess_orb()

subroutine thesky_coordinates::precess_orb ( real(double), intent(in) jd1,
real(double), intent(in) jd2,
real(double), intent(inout) i,
real(double), intent(inout) o1,
real(double), intent(inout) o2 )

Compute the precession of the equinoxes in orbital elements.

Parameters
jd1Original Julian day
jd2Target Julian day
iInclination (output)
o1Argument of perihelion (output)
o2Longitude of ascending node (output)
See also
Meeus, Astronomical Algorithms, 1998, Ch.24

Definition at line 341 of file coordinates.f90.

Referenced by thesky_asteroids::asteroid_elements().

◆ precess_xyz()

subroutine thesky_coordinates::precess_xyz ( real(double), intent(in) jd1,
real(double), intent(in) jd2,
real(double), intent(inout) x,
real(double), intent(inout) y,
real(double), intent(inout) z )

Compute the precession of the equinoxes in rectangular coordinates, from jd1 to jd2.

Parameters
jd1Original Julian day
jd2Target Julian day
xGeocentric rectangular X (output)
yGeocentric rectangular Y (output)
zGeocentric rectangular Z (output)
See also
Meeus, Astronomical Algorithms, 1998, Ch. 21

Definition at line 148 of file coordinates.f90.

Referenced by thesky_comets::cometgc().

◆ rect_2_spher()

subroutine thesky_coordinates::rect_2_spher ( real(double), intent(in) x,
real(double), intent(in) y,
real(double), intent(in) z,
real(double), intent(out) l,
real(double), intent(out) b,
real(double), intent(out) r )

Convert rectangular coordinates x,y,z to spherical coordinates l,b,r.

Parameters
xRectangular x coordinate (same unit as r)
yRectangular y coordinate (same unit as r)
zRectangular z coordinate (same unit as r)
lLongitude (rad) (output)
bLatitude (rad) (output)
rDistance (same unit as x,y,z) (output)

Definition at line 388 of file coordinates.f90.

Referenced by thesky_planets::planet_position(), and thesky_planets::saturnphys().

◆ refi()

real(double) function thesky_coordinates::refi ( real(double), intent(in) r,
real(double), intent(in) n,
real(double), intent(in) dndr )

The refraction integrand for atmospheric_refraction_apparent()

Parameters
rThe current distance from the centre of the Earth in metres
nThe refractive index at R
dndrThe rate the refractive index is changing at R
Return values
refiThe integrand of the refraction function

Definition at line 1302 of file coordinates.f90.

References refi().

Referenced by atmospheric_refraction_apparent(), and refi().

Here is the call graph for this function:

◆ refract()

real(double) function thesky_coordinates::refract ( real(double), intent(in) alt,
real(double), intent(in), optional press,
real(double), intent(in), optional temp )

Compute the atmospheric refraction for a given true altitude. You should add the result to the uncorrected altitude in order to obtain the observed altitude. Return 0 if alt + refract < -0.3°.

Parameters
altTrue (computed) altitude (rad)
pressAir pressure (hPa; optional)
tempAir temperature (degrees Celsius; optional)
Return values
refractRefraction in altitude (rad). You should add the result to the uncorrected altitude.
See also
  • T. Saemundsson, "Atmospheric refraction", Sky & Telescope vol.72, p.70 (1986), converted to radians;
  • Meeus (1998), Eq. 16.4 ff.

Definition at line 910 of file coordinates.f90.

References refract().

Referenced by thesky_planets::planet_position(), thesky_planets::planet_position_la(), refract(), refract_sp(), thesky_riset::riset(), thesky_riset::riset_ipol(), and thesky_sun::sunpos_la().

Here is the call graph for this function:

◆ refract_sp()

real function thesky_coordinates::refract_sp ( real, intent(in) alt,
real, intent(in), optional press,
real, intent(in), optional temp )

Compute the atmospheric refraction for a given true altitude, using single-precision values. This is a wrapper for refract().

Parameters
altTrue (computed) altitude (rad)
pressAir pressure (hPa; optional)
tempAir temperature (degrees Celsius; optional)
Return values
refractRefraction in altitude (rad). You should add the result to the uncorrected altitude.
See also
refract() for details.

Definition at line 951 of file coordinates.f90.

References refract(), and refract_sp().

Referenced by refract_sp().

Here is the call graph for this function:

◆ stratosphere_model()

subroutine thesky_coordinates::stratosphere_model ( real(double), intent(in) rt,
real(double), intent(in) tt,
real(double), intent(in) nt,
real(double), intent(in) a,
real(double), intent(in) r,
real(double), intent(out) n,
real(double), intent(out) dndr )

Stratosphere model for atmospheric_refraction_apparent()

Parameters
rtThe height of the tropopause from the centre of the Earth in metres
ttThe temperature at the tropopause in Kelvin
ntThe refractive index at the tropopause
aConstant of the atmospheric model = G*MD/R
rThe current distance from the centre of the Earth in metres
nThe refractive index at R (output)
dndrThe rate the refractive index is changing at R (output)

Definition at line 1356 of file coordinates.f90.

Referenced by atmospheric_refraction_apparent().

◆ troposphere_model()

subroutine thesky_coordinates::troposphere_model ( real(double), intent(in) r0,
real(double), intent(in) t0,
real(double), dimension(10), intent(in) a,
real(double), intent(in) r,
real(double), intent(out) t,
real(double), intent(out) n,
real(double), intent(out) dndr )

Troposphere model for atmospheric_refraction_apparent()

Parameters
r0The height of the observer from the centre of the Earth
t0The temperature at the observer in Kelvin
aConstants defined at the observer
rThe current distance from the centre of the Earth in metres
tThe temperature at R in Kelvin (output)
nThe refractive index at R (output)
dndrThe rate the refractive index is changing at R (output)

Definition at line 1325 of file coordinates.f90.

Referenced by atmospheric_refraction_apparent().