HorizontalCoordinates

class sofia_redux.scan.coordinate_systems.horizontal_coordinates.HorizontalCoordinates(coordinates=None, unit='degree', copy=True)[source]

Bases: SphericalCoordinates

Initialize horizontal coordinates.

The horizontal coordinate system uses the observer’s local horizon as the fundamental plane defining the altitude (latitude) and longitude (azimuth) in spherical coordinates. The altitude (or elevation) is the angle between the object and the observer’s local horizon, usually measured between 0 and 90 degrees. The azimuth is the angle of the object around the horizon, measured from north to east.

Parameters:
coordinateslist or tuple or array-like or units.Quantity, optional

The coordinates used to populate the object during initialization. The first (0) value or index should represent longitudinal coordinates, and the second should represent latitude.

unitunits.Unit or str, optional

The angular unit for the spherical coordinates. The default is ‘degree’.

copybool, optional

Whether to explicitly perform a copy operation on the input coordinates when storing them into these coordinates. Note that it is extremely unlikely for the original coordinates to be passed in as a reference due to the significant checks performed on them.

Attributes Summary

az

Return the azimuth coordinates.

el

Return the elevation coordinates.

fits_latitude_stem

Return the FITS header latitude stem string.

fits_longitude_stem

Return the FITS header longitude stem string.

two_letter_code

Return the two-letter code representing the horizontal system.

za

Return the Zenith Angle (ZA).

Methods Summary

convert_horizontal_to_equatorial(horizontal, ...)

Convert horizontal coordinates to an equatorial frame.

copy()

Return a copy of the horizontal coordinates.

edit_header(header, key_stem[, alt])

Edit the header with horizontal coordinate information.

get_parallactic_angle(site)

Return the parallactic angle for a given site.

plot(*args, **kwargs)

Plot the coordinates.

set_az(azimuth)

Set the azimuth coordinates.

set_el(elevation)

Set the elevation coordinates.

set_za(zenith_angle)

Set the zenith angle.

setup_coordinate_system()

Setup the system for the coordinates.

to_equatorial(site, lst[, equatorial])

Return these coordinates as equatorial coordinates.

to_equatorial_offset(offset, position_angle)

Convert horizontal offsets to an equatorial offset.

Attributes Documentation

az

Return the azimuth coordinates.

Returns:
azimuthastropy.units.Quantity (float or numpy.ndarray)
el

Return the elevation coordinates.

Returns:
elevationastropy.units.Quantity (float or numpy.ndarray)
fits_latitude_stem

Return the FITS header latitude stem string.

Returns:
stemstr
fits_longitude_stem

Return the FITS header longitude stem string.

Returns:
stemstr
two_letter_code

Return the two-letter code representing the horizontal system.

Returns:
codestr
za

Return the Zenith Angle (ZA).

Returns:
astropy.units.Quantity (float or numpy.ndarray)

Methods Documentation

classmethod convert_horizontal_to_equatorial(horizontal, equatorial, site, lst)[source]

Convert horizontal coordinates to an equatorial frame.

Parameters:
horizontalHorizontalCoordinates
equatorialEquatorialCoordinates

The frame in which to hold the equatorial coordinates.

siteGeodeticCoordinates
lstastropy.units.Quantity (float or numpy.ndarray)
Returns:
None
copy()[source]

Return a copy of the horizontal coordinates.

Returns:
HorizontalCoordinates
edit_header(header, key_stem, alt='')[source]

Edit the header with horizontal coordinate information.

Parameters:
headerastropy.io.fits.header.Header

The header to modify.

key_stemstr

The name of the header key to update.

altstr, optional

The alternative coordinate system.

Returns:
None
get_parallactic_angle(site)[source]

Return the parallactic angle for a given site.

Parameters:
siteGeodeticCoordinates
Returns:
paastropy.units.Quantity (float or numpy.ndarray)
plot(*args, **kwargs)[source]

Plot the coordinates.

Parameters:
argsvalues

Optional positional parameters to pass into pyplot.plot.

kwargsdict, optional

Optional keyword arguments.

Returns:
None
set_az(azimuth)[source]

Set the azimuth coordinates.

Parameters:
azimuthastropy.units.Quantity (float or numpy.ndarray)
Returns:
None
set_el(elevation)[source]

Set the elevation coordinates.

Parameters:
elevationastropy.units.Quantity (float or numpy.ndarray)
Returns:
None
set_za(zenith_angle)[source]

Set the zenith angle.

Parameters:
zenith_angleastropy.units.Quantity (float or numpy.ndarray)
Returns:
None
setup_coordinate_system()[source]

Setup the system for the coordinates.

Returns:
None
to_equatorial(site, lst, equatorial=None)[source]

Return these coordinates as equatorial coordinates.

Parameters:
siteGeodeticCoordinates
lstastropy.units.Quantity
equatorialEquatorialCoordinates, optional
Returns:
EquatorialCoordinates
classmethod to_equatorial_offset(offset, position_angle, in_place=True)[source]

Convert horizontal offsets to an equatorial offset.

The offsets are updated in-place.

Parameters:
offsetCoordinate2D

The (x, y) horizontal offsets to rotate.

position_angleastropy.units.Quantity (float or numpy.ndarray)

The position angle as a scalar or array of shape (n,).

in_placebool, optional

If True, update offset in-place. Otherwise, return a new offset without modifying the original.

Returns:
equatorial_offsetCoordinate2D