Extension of astropy.coordinates


Overview

The grand package extends astropy.coordinates with geographic representations and local frames. It adds the ECEF and LTP frames as well as the GeodeticRepresentation and HorizontalRepresentation of coordinates.

Frames

The grand frames differ from base astropy.coordinates frames in several respects. First, they all inherit from a ExtendedCoordinateFrame instead of a BaseCoordinateFrame. Their coordinates data are read only. This prevents some severe bugs (e.g. see issue #9873). In addition, they support basic arithmetic operators (+, - and *).

Secondly, the grand frames are co-moving with the Earth. i.e. the obstime can be omitted, in which case it is inherited from the source frame during transforms. This differs from base astropy.coordinates frames where the observation time must be stated explicitly during a transform (E.g. see issue #8390).

If specified, the observation time must be an instance of an astropy Time object or any of its initializers, E.g. datetime, or str.

In addition, the frame classes defined in astropy.coordinates do not distinguish points from vectors when transforming between frames. The grand frames do. An automatic type inference occurs based on the data astropy.units, as following.

Warning

If the frame data are homogeneous to a length, there are assumed to transform as a point. Otherwise, they transform as a vector quantity.

ECEF frame

The ECEF class (Earth-Centered Earth-Fixed) is a wrapper for geocentric frames. It behaves as the ITRS frame but co-moving with the Earth, i.e. the obstime can be omitted.

class grand.ECEF(*args, obstime: Union[datetime.datetime, astropy.time.core.Time, str, None] = None, **kwargs)

Earth-Centered Earth-Fixed frame, co-moving with the Earth.

Note

By default, a CartesianRepresentation is expected for the coordinates data, i.e. x, y and z. For example, the following represents a point located at the origin.

>>> origin = ECEF(x = 0 * u.m, y = 0 * u.m, z = 0 * u.m)
property earth_location

The data in this frame as an EarthLocation.

property obstime

The observation time. Defaults to None.

LTP frame

The LTP class (Local Tangent Plane coordinates) allows to define local frames tangent to the WGS84 ellipsoid. The orientation property can be used to define the frame axes, along the directions of the local magnetic field. E.g. the following defines a local magnetic North, East, Down (NED) frame centered on Greenwich, as:

>>> from astropy.coordinates import EarthLocation
>>> ltp = LTP(location=EarthLocation.of_site('greenwich'),
...           orientation='NED')

By default magnetic North, West, Up (NWU) coordinates are used. If not explicitly specified the magnetic declination is computed using the built-in geomagnet module. Alternatively, geographic coordinates can be used as well by setting magnetic to False.

Note

One must always specify the local frame origin, as the location parameter.

class grand.LTP(*args, location: Union[EarthLocation, ECEF, LTP, None] = None, orientation: Optional[Sequence[str]] = None, magnetic: Optional[bool] = None, declination: Optional[astropy.units.quantity.Quantity] = None, rotation: Optional[grand.tools.coordinates.transform.Rotation] = None, obstime: Union[datetime, Time, str, None] = None, **kwargs)

Local frame tangent to the WGS84 ellipsoid and oriented along cardinal directions.

Note

By default, a CartesianRepresentation is expected for the coordinates data, i.e. x, y and z.

rotated(rotation: grand.tools.coordinates.transform.Rotation, copy: bool = True) → LTP

Get a rotated version of this frame.

property earth_location

The data in this frame as an EarthLocation.

property magnetic

Use the magnetic north instead of the geographic one (default: True).

property orientation

The cardinal directions of the x, y, and z axis (default: ENU).

property obstime

The observation time.

Extended frame

The ExtendedCoordinateFrame is an extension of the astropy BaseCoordinateFrame but with read only data. This prevents some type of severe bugs that could occur with astropy frames when data are modified in-place (E.g. see issue #9873).

Furthermore the ExtendedCoordinateFrame also implements the addition (+), subtraction (-) and multiplication (*) arithmetic operators. Those are forwarded to the frame coordinates.

Note

When adding or subtracting two frames the resulting system is the one of the left hand side. If needed, the coordinates of the right hand side are transformed to the system of the left hand side before being added or subtracted. Consequently, the addition of two frames is not strictly commutative. Note that one can also add or subtract a frame with a BaseRepresentation in which case the representation is assumed to be in the frame system.

class grand.ExtendedCoordinateFrame(*args, **kwargs)

A coordinates frame with immutable data supporting extra arithmetic operators.

This class inherits from a BaseCoordinateFrame. Most functionalities are identical and are not documented here.

Warning

The data managed by the frame and by any of its representations are enforced to be read only in order to prevent in-place modifications. However, contrary to a astropy.coordinates.BaseCoordinateFrame it is allowed to (re)set the data to another representation, e.g. as:

>>> frame.data = CartesianRepresentation(0, 0, 0)

Representations

The grand package explictly adds two geographic representations based on the WGS84 ellipsoid. Note that although astropy uses the corresponding transforms it does not explictly bundle them as coordinates representations.

Geodetic coordinates

A GeodeticRepresentation is analog to a SphericalRepresentation but mapping to the WGS84 ellipsoid instead of a sphere. It allows to represent ECEF coordinates using a geodetic datum instead of the default astropy.coordinates.Cartesian, E.g. as:

>>> point = ECEF(latitude=45 * u.deg, longitude=3 * u.deg, height=0.5 * u.m,
...              representation_type='geodetic')
class grand.GeodeticRepresentation(latitude: astropy.coordinates.angles.Latitude, longitude: astropy.coordinates.angles.Longitude, height: astropy.units.quantity.Quantity = <Quantity 0. m>, copy: bool = True)

Geodetic coordinates representation w.r.t. the WGS84 ellipsoid.

Note

The latitude angle is measured clockwise, w.r.t. the xOy plane. The longitude angle is measured counter-clockwise, w.r.t. the x-axis. The height is w.r.t. the WGS84 ellipsoid. It defaults to 0 meters if ommitted.

classmethod from_cartesian(cart: astropy.coordinates.representation.CartesianRepresentation) → grand.tools.coordinates.representation.GeodeticRepresentation

Generate a Geodetic representation from a Cartesian one.

to_cartesian() → astropy.coordinates.representation.CartesianRepresentation

Generate a Cartesian representation from a Geodetic one.

Horizontal coordinates

The HorizontalRepresentation allows to represent a unit vector, e.g. a direction, in a local LTP frame using a horizontal coordinates system. For example, the following defines a unit vector pointing upwards at Greenwich:

>>> from astropy.coordinates import EarthLocation
>>> up = LTP(location=EarthLocation.of_site('greenwich'),
...          representation_type='horizontal',
...          azimuth = 0 * u.deg, elevation = 90 * u.deg)
class grand.HorizontalRepresentation(azimuth: Union[astropy.coordinates.angles.Longitude, astropy.units.quantity.Quantity, str], elevation: Union[astropy.coordinates.angles.Latitude, astropy.units.quantity.Quantity, str], copy: bool = True)

Horizontal angular representation, for unit vectors.

Note

The azimuth angle is measured clockwise, w.r.t. the y axis. The elevation angle is measured clockwise, w.r.t. the xOy plane.

classmethod from_cartesian(cart: astropy.coordinates.representation.CartesianRepresentation) → grand.tools.coordinates.representation.HorizontalRepresentation

Generate a Horizontal angular representation from a Cartesian unit vector.

Warning

The Cartesian unit vector must be dimensioneless. Though it is not checked if the norm of the vector is indeed unity.

to_cartesian() → astropy.coordinates.representation.CartesianRepresentation

Generate a Cartesian unit vector from this Horizontal angular representation.