SpectralCoord

class astropy.coordinates.SpectralCoord(value, unit=None, observer=None, target=None, radial_velocity=None, redshift=None, **kwargs)[source]

Bases: SpectralQuantity

A spectral coordinate with its corresponding unit.

Note

The SpectralCoord class is new in Astropy v4.1 and should be considered experimental at this time. Note that we do not fully support cases where the observer and target are moving relativistically relative to each other, so care should be taken in those cases. It is possible that there will be API changes in future versions of Astropy based on user feedback. If you have specific ideas for how it might be improved, please let us know on the astropy-dev mailing list or at http://feedback.astropy.org.

Parameters:
valuendarray or Quantity or SpectralCoord

Spectral values, which should be either wavelength, frequency, energy, wavenumber, or velocity values.

unitastropy:unit-like

Unit for the given spectral values.

observerBaseCoordinateFrame or SkyCoord, optional

The coordinate (position and velocity) of observer. If no velocities are present on this object, the observer is assumed to be stationary relative to the frame origin.

targetBaseCoordinateFrame or SkyCoord, optional

The coordinate (position and velocity) of target. If no velocities are present on this object, the target is assumed to be stationary relative to the frame origin.

radial_velocityQuantity [:ref: ‘speed’], optional

The radial velocity of the target with respect to the observer. This can only be specified if redshift is not specified.

redshiftpython:float, optional

The relativistic redshift of the target with respect to the observer. This can only be specified if radial_velocity cannot be specified.

doppler_restQuantity, optional

The rest value to use when expressing the spectral value as a velocity.

doppler_conventionpython:str, optional

The Doppler convention to use when expressing the spectral value as a velocity.

Attributes Summary

observer

The coordinates of the observer.

quantity

Convert the SpectralCoord to a Quantity.

radial_velocity

Radial velocity of target relative to the observer.

redshift

Redshift of target relative to observer.

target

The coordinates of the target being observed.

Methods Summary

replicate([value, unit, observer, target, ...])

Return a replica of the SpectralCoord, optionally changing the values or attributes.

to_rest()

Transforms the spectral axis to the rest frame.

with_observer_stationary_relative_to(frame)

A new SpectralCoord with the velocity of the observer altered, but not the position.

with_radial_velocity_shift([target_shift, ...])

Apply a velocity shift to this spectral coordinate.

Attributes Documentation

observer

The coordinates of the observer.

If set, and a target is set as well, this will override any explicit radial velocity passed in.

Returns:
BaseCoordinateFrame

The astropy coordinate frame representing the observation.

quantity

Convert the SpectralCoord to a Quantity. Equivalent to self.view(u.Quantity).

Returns:
Quantity

This object viewed as a Quantity.

radial_velocity

Radial velocity of target relative to the observer.

Returns:
Quantity [:ref: ‘speed’]

Radial velocity of target.

Notes

This is different from the .radial_velocity property of a coordinate frame in that this calculates the radial velocity with respect to the observer, not the origin of the frame.

redshift

Redshift of target relative to observer. Calculated from the radial velocity.

Returns:
astropy.units.Quantity

Redshift of target.

target

The coordinates of the target being observed.

If set, and an observer is set as well, this will override any explicit radial velocity passed in.

Returns:
BaseCoordinateFrame

The astropy coordinate frame representing the target.

Methods Documentation

replicate(value=None, unit=None, observer=None, target=None, radial_velocity=None, redshift=None, doppler_convention=None, doppler_rest=None, copy=False)[source]

Return a replica of the SpectralCoord, optionally changing the values or attributes.

Note that no conversion is carried out by this method - this keeps all the values and attributes the same, except for the ones explicitly passed to this method which are changed.

If copy is set to True then a full copy of the internal arrays will be made. By default the replica will use a reference to the original arrays when possible to save memory.

Parameters:
valuendarray or Quantity or SpectralCoord, optional

Spectral values, which should be either wavelength, frequency, energy, wavenumber, or velocity values.

unitastropy:unit-like

Unit for the given spectral values.

observerBaseCoordinateFrame or SkyCoord, optional

The coordinate (position and velocity) of observer.

targetBaseCoordinateFrame or SkyCoord, optional

The coordinate (position and velocity) of target.

radial_velocityQuantity [:ref: ‘speed’], optional

The radial velocity of the target with respect to the observer.

redshiftpython:float, optional

The relativistic redshift of the target with respect to the observer.

doppler_restQuantity, optional

The rest value to use when expressing the spectral value as a velocity.

doppler_conventionpython:str, optional

The Doppler convention to use when expressing the spectral value as a velocity.

copybool, optional

If True, and value is not specified, the values are copied to the new SkyCoord - otherwise a reference to the same values is used.

Returns:
scSpectralCoord object

Replica of this object

to_rest()[source]

Transforms the spectral axis to the rest frame.

with_observer_stationary_relative_to(frame, velocity=None, preserve_observer_frame=False)[source]

A new SpectralCoord with the velocity of the observer altered, but not the position.

If a coordinate frame is specified, the observer velocities will be modified to be stationary in the specified frame. If a coordinate instance is specified, optionally with non-zero velocities, the observer velocities will be updated so that the observer is co-moving with the specified coordinates.

Parameters:
framepython:str, BaseCoordinateFrame or SkyCoord

The observation frame in which the observer will be stationary. This can be the name of a frame (e.g. ‘icrs’), a frame class, frame instance with no data, or instance with data. This can optionally include velocities.

velocityQuantity or CartesianDifferential, optional

If frame does not contain velocities, these can be specified as a 3-element Quantity. In the case where this is also not specified, the velocities default to zero.

preserve_observer_framebool

If True, the final observer frame class will be the same as the original one, and if False it will be the frame of the velocity reference class.

Returns:
new_coordSpectralCoord

The new coordinate object representing the spectral data transformed based on the observer’s new velocity frame.

with_radial_velocity_shift(target_shift=None, observer_shift=None)[source]

Apply a velocity shift to this spectral coordinate.

The shift can be provided as a redshift (float value) or radial velocity (Quantity with physical type of ‘speed’).

Parameters:
target_shiftpython:float or Quantity [:ref: ‘speed’]

Shift value to apply to current target.

observer_shiftpython:float or Quantity [:ref: ‘speed’]

Shift value to apply to current observer.

Returns:
SpectralCoord

New spectral coordinate with the target/observer velocity changed to incorporate the shift. This is always a new object even if target_shift and observer_shift are both None.