Horizontal coordinates

Convert hour angle and declination into horizontal coordinates

PyAstronomy.pyasl.hadec2altaz(ha, dec, lat, ws=False, radian=False)

Convert hour angle and declination into horizon (alt/az) coordinates.

Parameters
hafloat or array

Local apparent hour angle in DEGREES.

decfloat or array

Local apparent declination in DEGREES.

latfloat or array

Local latitude in DEGREES.

radianboolean, optional

If True, the result is returned in radian instead of in degrees (default is False).

wsboolean, optional

Set this to True, if the azimuth shall be measured West from South. Default is to measure azimuth East from North.

Returns
Altitudelist

A list holding the Local Apparent Altitude [deg].

Apparent Azimuthlist

The Local Apparent Azimuth [deg].

Notes

Note

This function was ported from the IDL Astronomy User’s Library.

IDL - Documentation
NAME:

HADEC2ALTAZ

PURPOSE:

Converts Hour Angle and Declination to Horizon (alt-az) coordinates.

EXPLANATION:

Can deal with NCP/SCP singularity. Intended mainly to be used by program EQ2HOR

CALLING SEQUENCE:

HADEC2ALTAZ, ha, dec, lat ,alt ,az [ /WS ]

INPUTS

ha - the local apparent hour angle, in DEGREES, scalar or vector dec - the local apparent declination, in DEGREES, scalar or vector lat - the local latitude, in DEGREES, scalar or vector

OUTPUTS

alt - the local apparent altitude, in DEGREES. az - the local apparent azimuth, in DEGREES, all results in double

precision

OPTIONAL KEYWORD INPUT:
/WS - Set this keyword for the output azimuth to be measured West from

South. The default is to measure azimuth East from North.

EXAMPLE:

What were the apparent altitude and azimuth of the sun when it transited the local meridian at Pine Bluff Observatory (Lat=+43.07833 degrees) on April 21, 2002? An object transits the local meridian at 0 hour angle. Assume this will happen at roughly 1 PM local time (18:00 UTC).

IDL> jdcnv, 2002, 4, 21, 18., jd ; get rough Julian date to determine

;Sun ra, dec.

IDL> sunpos, jd, ra, dec IDL> hadec2altaz, 0., dec, 43.078333, alt, az

===> Altitude alt = 58.90

Azimuth az = 180.0

REVISION HISTORY:

Written Chris O’Dell Univ. of Wisconsin-Madison May 2002

Example:

from __future__ import print_function, division
from PyAstronomy import pyasl
import datetime
import numpy as np

# Hour angle 0. means transiting the local meridian.
ha = 0.
# Declination of object
dec = 30.
# Latitude of the observer (here Hamburger Sternwarte)
lat = +53.48
print("Get altitude and azimuth of object in DEGREES")
print(pyasl.hadec2altaz(ha, dec, lat))

# List of coordinates
ha = np.arange(0., 20., 5.)
dec = np.arange(30., 50., 5.)
lat = np.zeros(dec.size)+53.48
print()
print("Get altitude and azimuth for a list of objects from same observer latitude")
altaz = pyasl.hadec2altaz(ha, dec, lat)
print("alt: ", altaz[0])
print("az: ", altaz[1])

Convert celestial coordinates (RA/DEC) into local horizon coordinates (ALT/AZ)

PyAstronomy.pyasl.eq2hor(jd, ra, dec, observatory=None, lon=None, lat=None, alt=None, B1950=False, precess=True, nutate=True, aberration=True, refract=True)

Convert celestial coordinates (RA/DEC) to local horizon coordinates (ALT/AZ).

This routine is typically accurate to about 1 arcsec. It considers Earth’s precession, nutation, aberration, and refraction (if keywords are True).

Parameters
jdfloat or array

The Julian date(s)

rafloat or array

The right ascension in DEGREES.

decfloat or array

The declination in DEGREES.

observatorystring, {HS}, optional

A string identifying the observatory. If given, the observer’s longitude, latitude, and altitude are set automatically (and must not be given separately then).

lonfloat

East longitude of the observer in DEGREES. Specify West longitude with a negative sign. Default is the longitude of the Hamburger Sternwarte.

latfloat

Latitude of the observer in DEGREES. Default is the latitude of the Hamburger Sternwarte.

altfloat

Altitude of the observer in METER. Default is the altitude of the Hamburger Sternwarte.

B1950boolean, optional

If True, your RA and DEC coordinates are given for epoch B1950 FK4. If False, RA and DEC are given for epoch J2000 FK5. Default is FALSE.

precessboolean, optional

If True (default), Earth’s precess motion is considered in the calculations.

nutateboolean, optional

If True (default), Earth’s nutation is considered in the calculations.

aberrationboolean, optional

If True (default), the annual aberration is considered in the calculations.

refractionboolean, optional

If True, the atmospheric refraction is considered in the calculations.

Returns
Altitudefloat or array

The altitude in degrees.

Azimuthfloat or array

The azimuth in degrees (measured East from North).

Hour anglefloat or array

The hour angle in degrees.

Notes

Note

This function was ported from the IDL Astronomy User’s Library.

IDL - Documentation

NAME:

EQ2HOR

PURPOSE:

Convert celestial (ra-dec) coords to local horizon coords (alt-az).

CALLING SEQUENCE:

eq2hor, ra, dec, jd, alt, az, [ha, LAT= , LON= , /WS, OBSNAME= , $

/B1950 , PRECESS_= 0, NUTATE_= 0, REFRACT_= 0, $ ABERRATION_= 0, ALTITUDE= , /VERBOSE, _EXTRA= ]

DESCRIPTION:

This code calculates horizon (alt,az) coordinates from equatorial (ra,dec) coords. It is typically accurate to about 1 arcsecond or better (I have checked the output against the publicly available XEPHEM software). It performs precession, nutation, aberration, and refraction corrections. The perhaps best thing about it is that it can take arrays as inputs, in all variables and keywords EXCEPT Lat, lon, and Altitude (the code assumes these aren’t changing), and uses vector arithmetic in every calculation except when calculating the precession matrices.

INPUT VARIABLES:
RARight Ascension of object (J2000) in degrees (FK5); scalar or

vector.

Dec : Declination of object (J2000) in degrees (FK5), scalar or vector. JD : Julian Date [scalar or vector]

Note: if RA and DEC are arrays, then alt and az will also be arrays.

If RA and DEC are arrays, JD may be a scalar OR an array of the same dimensionality.

OPTIONAL INPUT KEYWORDS:

lat : north geodetic latitude of location in degrees lon : EAST longitude of location in degrees (Specify west longitude

with a negative sign.)

WSSet this to get the azimuth measured westward from south (not

East of North).

obsname: Set this to a valid observatory name to be used by the

astrolib OBSERVATORY procedure, which will return the latitude and longitude to be used by this program.

B1950Set this if your ra and dec are specified in B1950, FK4

coordinates (instead of J2000, FK5)

precessSet this to 1 to force precession [default], 0 for no

precession correction

nutate : Set this to 1 to force nutation [default], 0 for no nutation. aberration : Set this to 1 to force aberration correction [default],

0 for no correction.

refractSet to 1 to force refraction correction [default], 0 for no

correction.

altitude: The altitude of the observing location, in meters. [default=0]. verbose: Set this for verbose output. The default is verbose=0. extra: This is for setting TEMPERATURE or PRESSURE explicitly, which are

used by CO_REFRACT to calculate the refraction effect of the atmosphere. If you don’t set these, the program will make an intelligent guess as to what they are (taking into account your altitude). See CO_REFRACT for more details.

OUTPUT VARIABLES: (all double precision)

alt : altitude (in degrees) az : azimuth angle (in degrees, measured EAST from NORTH, but see

keyword WS above.)

ha : hour angle (in degrees) (optional)

DEPENDENCIES:

NUTATE, PRECESS, OBSERVATORY, SUNPOS, ADSTRING() CO_NUTATE, CO_ABERRATION, CO_REFRACT, ALTAZ2HADEC, SETDEFAULTVALUE

BASIC STEPS

Apply refraction correction to find apparent Alt. Calculate Local Mean Sidereal Time Calculate Local Apparent Sidereal Time Do Spherical Trig to find apparent hour angle, declination. Calculate Right Ascension from hour angle and local sidereal time. Nutation Correction to Ra-Dec Aberration correction to Ra-Dec

Precess Ra-Dec to current equinox.

ORRECTIONS I DO NOT MAKE:
  • Deflection of Light by the sun due to GR. (typically milliarcseconds,

    can be arseconds within one degree of the sun)

  • The Effect of Annual Parallax (typically < 1 arcsecond)

  • and more (see below)

TO DO
  • Better Refraction Correction. Need to put in wavelength dependence,

and integrate through the atmosphere.
  • Topocentric Parallax Correction (will take into account elevation of the observatory)

  • Proper Motion (but this will require crazy lookup tables or something).
    • Difference between UTC and UT1 in determining LAST – is this important?

    • Effect of Annual Parallax (is this the same as topocentric Parallax?)

  • Polar Motion
    • Better connection to Julian Date Calculator.

EXAMPLE

Find the position of the open cluster NGC 2264 at the Effelsburg Radio Telescope in Germany, on June 11, 2023, at local time 22:00 (METDST). The inputs will then be:

Julian Date = 2460107.250 Latitude = 50d 31m 36s Longitude = 06h 51m 18s Altitude = 369 meters RA (J2000) = 06h 40m 58.2s Dec(J2000) = 09d 53m 44.0s

IDL> eq2hor, ten(6,40,58.2)*15., ten(9,53,44), 2460107.250d, alt, az, $
lat=ten(50,31,36), lon=ten(6,51,18), altitude=369.0, /verb, $

pres=980.0, temp=283.0

The program produces this output (because the VERBOSE keyword was set)

Latitude = +50 31 36.0 Longitude = +06 51 18.0 Julian Date = 2460107.250000 Ra, Dec: 06 40 58.2 +09 53 44.0 (J2000) Ra, Dec: 06 42 15.7 +09 52 19.2 (J2023.4422) Ra, Dec: 06 42 13.8 +09 52 26.9 (fully corrected) LMST = +11 46 42.0 LAST = +11 46 41.4 Hour Angle = +05 04 27.6 (hh:mm:ss) Az, El = 17 42 25.6 +16 25 10.3 (Apparent Coords) Az, El = 17 42 25.6 +16 28 22.8 (Observer Coords)

Compare this with the result from XEPHEM: Az, El = 17h 42m 25.6s +16d 28m 21s

This 1.8 arcsecond discrepancy in elevation arises primarily from slight differences in the way I calculate the refraction correction from XEPHEM, and is pretty typical.

AUTHOR:
Chris O’Dell

Univ. of Wisconsin-Madison

Observational Cosmology Laboratory Email: odell@cmb.physics.wisc.edu

Revision History:

August 2012 Use Strict_Extra to flag spurious keywords W. Landsman

Example:

from __future__ import print_function, division
from PyAstronomy import pyasl
import matplotlib.pylab as plt
import datetime
import numpy as np

# Convert calendar date to JD
# use the datetime package
jd = datetime.datetime(2013, 4, 16)
jd = pyasl.jdcnv(jd)
# Specific RA and DEC
ra = 10.
dec = 30.
print()
print("Get horizontal coordinates (alt, az, ha) from JD, RA,")
print("  and DEC for the Hamburger Sternwarte")
print(pyasl.eq2hor(jd, ra, dec, observatory="HS"))

print()
print("From a list of Julian dates ...")
jds = np.arange(jd, jd+1, .2)
ras = np.zeros(jds.size) + ra
decs = np.zeros(jds.size) + dec
alt, az, ha = pyasl.eq2hor(jds, ras, decs, lon=-70.4042, lat=-24.6272, alt=2635.)

for i in range(alt.size):
    print("JD = %g : alt = % g,  az = % g,  ha = % g" % (jds[i], alt[i], az[i], ha[i]))


print()
print("For one object and different times at the VLT...")
jds = np.arange(jd-.25, jd+.25, .01)
ras = np.zeros(jds.size) + 130.
decs = np.zeros(jds.size) - 30.
res = pyasl.eq2hor(jds, ras, decs, lon=-70.4042, lat=-24.6272, alt=2635.)

plt.plot(jds, res[0])
plt.xlabel("Julian date")
plt.ylabel("Altitude [deg]")
plt.show()