Package pywcs :: Class WCSBase
[frames] | no frames]

Class WCSBase

source code

object --+
         |
        WCSBase
Known Subclasses:

WCS(header, key=' ', relax=False)

WCS objects convert between pixel and world coordinates, based on the WCS settings in a FITS file.

Instance Methods
 
__copy__()
Creates a deep copy of the WCS object.
source code
 
__init__(header, key=' ', relax=False)
x.__init__(...) initializes x; see x.__class__.__doc__ for signature
source code
a new object with type S, a subtype of T
__new__(T, S, ...) source code
int
celfix()
Translates AIPS-convention celestial projection types, -NCP and -GLS.
source code
 
copy()
Creates a deep copy of the WCS object.
source code
int
cylfix()
Fixes WCS keyvalues for malformed cylindrical projections.
source code
int
datfix()
Translates the old DATE-OBS date format to year-2000 standard form (yyyy-mm-ddThh:mm:ss) and derives MJD-OBS from it if not already set.
source code
dict
fix(translate_units='', naxis=0)
Applies all of the corrections handled separately by datfix, unitfix, celfix, spcfix and cylfix.
source code
bool
has_cdi_ja()
Returns True if CDi_ja is present.
source code
bool
has_crotaia()
Returns True if CROTAia is present.
source code
bool
has_pci_ja()
Returns True if PCi_ja is present.
source code
dict
mix(mixpix, mixcel, vspan, vstep, viter, world, pixcrd)
Given either the celestial longitude or latitude plus an element of the pixel coordinate, solves for the remaining elements by iterating on the unknown celestial coordinate element using s2p.
source code
dict
p2s(pixcrd)
Converts pixel to world coordinates.
source code
 
print_contents()
Print the contents of the WCS object to stdout.
source code
dict
s2p(world)
Transforms world coordinates to pixel coordinates.
source code
 
set()
Sets up a WCS object for use according to information supplied within it.
source code
int
spcfix()
Translates AIPS-convention spectral coordinate types.
source code
 
sptr(ctype, i=-1)
Translates the spectral axis in a WCS object.
source code
int
unitfix(translate_units='')
Translates non-standard CUNITia keyvalues.
source code

Inherited from object: __delattr__, __getattribute__, __hash__, __reduce__, __reduce_ex__, __repr__, __setattr__, __str__

Properties
string alt
Character code for alternate coordinate descriptions.
array[2][2] of double cd
CDi_ja linear transformation matrix.
array[naxis] of double cdelt
Coordinate increments (CDELTia) for each coord axis.
list of strings cname
A list of the coordinate axis names, CNAMEia.
array[naxis] of int colax
An array recording the column numbers for each axis in a pixel list.
  colnum
Where the coordinate representation is associated with an image-array column in a FITS binary table, this variable may be used to record the relevant column number.
array[naxis] of double crder
The random error in each coordinate axis, CRDERia.
array[2][2] of double crota
CROTAia keyvalues for each coordinate axis.
array[naxis] of double crpix
Coordinate reference pixels (CRPIXja) for each pixel axis.
array[naxis] of double crval
Coordinate reference values (CRVALia) for each coordinate axis.
array[naxis] of double csyer
The systematic error in the coordinate value axes, CSYERia.
list of strings ctype
List of CTYPEia keyvalues.
int cubeface
Index into the pixcrd (pixel coordinate) array for the CUBEFACE axis.
list of strings cunit
List of CUNITia keyvalues which define the units of measurement of the CRVALia, CDELTia and CDi_ja keywords.
string dateavg
Representative mid-point of the date of observation in ISO format, yyyy-mm-ddThh:mm:ss.
string dateobs
Start of the date of observation in ISO format, yyyy-mm-ddThh:mm:ss.
float equinox
The equinox associated with dynamical equatorial or ecliptic coordinate systems, EQUINOXa (or EPOCH in older headers).
int lat
The index into the world coordinate array containing latitude values.
float latpole
The native latitude of the celestial pole, LATPOLEa (deg).
int lng
The index into the world coordinate array containing longitude values.
float lonpole
The native longitude of the celestial pole, LONPOLEa (deg).
float mjdavg
Modified Julian Date (MJD = JD - 2400000.5), MJD-AVG, corresponding to DATE-AVG.
float mjdobs
Modified Julian Date (MJD = JD - 2400000.5), MJD-OBS, corresponding to DATE-OBS.
string name
The name given to the coordinate representation WCSNAMEa.
int naxis
The number of axes (pixel and coordinate), given by the NAXIS or WCSAXESa keyvalues.
array[3] of double obsgeo
Location of the observer in a standard terrestrial reference frame, OBSGEO-X, OBSGEO-Y, OBSGEO-Z (in metres).
array[2][2] of double pc
The PCi_ja (pixel coordinate) transformation matrix.
list of tuples ps
PSi_ma keywords for each i and m.
list of tuples pv
PVi_ma keywords for each i and m.
string radesys
The equatorial or ecliptic coordinate system type, RADESYSa.
float restfrq
Rest frequency (Hz) from RESTFRQa.
float restwav
Rest wavelength (m) from RESTWAVa.
int spec
The index containing the spectral axis values.
string specsys
Spectral reference frame (standard of rest), SPECSYSa
string ssysobs
The actual spectral reference frame in which there is no differential variation in the spectral coordinate across the field-of-view, SSYSOBSa.
string ssyssrc
The spectral reference frame (standard of rest) in which the redshift was measured, SSYSSRCa.
float velangl
The angle in degrees that should be used to decompose an observed velocity into radial and transverse components.
float velosys
The relative radial velocity (m/s) between the observer and the selected standard of rest in the direction of the celestial reference coordinate, VELOSYSa.
float zsource
The redshift, ZSOURCEa, of the source.

Inherited from object: __class__

Method Details

__init__(header, key=' ', relax=False)
(Constructor)

source code 

x.__init__(...) initializes x; see x.__class__.__doc__ for signature

Parameters:
  • header (PyFITS header object or string) - A PyFITS header object or a string containing the raw FITS header data.
  • key (string) - The key referring to a particular WCS transform in the header. This may be either ' ' or 'A'-'Z' and corresponds to the "a" part of "CTYPEia".
  • relax (bool) - Degree of permissiveness:
    • False: Recognize only FITS keywords defined by the published WCS standard.
    • True: Admit all recognized informal extensions of the WCS standard.
Raises:
  • MemoryError - Memory allocation failed.
  • ValueError - Invalid key.
  • KeyError - Key not found in FITS header.
Overrides: object.__init__

__new__(T, S, ...)

source code 
Returns: a new object with type S, a subtype of T
Overrides: object.__new__

celfix()

source code 

Translates AIPS-convention celestial projection types, -NCP and -GLS.

Returns: int
0 for success; -1 if no change required.

cylfix()

source code 

Fixes WCS keyvalues for malformed cylindrical projections.

Returns: int
0 for success; -1 if no change required.

datfix()

source code 

Translates the old DATE-OBS date format to year-2000 standard form (yyyy-mm-ddThh:mm:ss) and derives MJD-OBS from it if not already set. Alternatively, if mjdobs is set and dateobs isn't, then datfix derives dateobs from it. If both are set but disagree by more than half a day then ValueError is raised.

Returns: int
0 for success; -1 if no change required.

fix(translate_units='', naxis=0)

source code 

Applies all of the corrections handled separately by datfix, unitfix, celfix, spcfix and cylfix.

Parameters:
  • translate_units (string) - Do potentially unsafe translations of non-standard unit strings.

    Although "S" is commonly used to represent seconds, its translation to "s" is potentially unsafe since the standard recognizes "S" formally as Siemens, however rarely that may be used. The same applies to "H" for hours (Henry), and "D" for days (Debye).

    This string controls what to do in such cases, and is case-insensitive.

    • If the string contains "s", translate "S" to "s".
    • If the string contains "h", translate "H" to "h".
    • If the string contains "d", translate "D" to "d".

    Thus '' doesn't do any unsafe translations, whereas 'shd' does all of them.

  • naxis (array[naxis] of int) - Image axis lengths. If this array pointer is set to zero, then cylfix will not be invoked.
Returns: dict
A dictionary containing the following keys, each referring to a status string for each of the sub-fix functions that were called: datfix, unitfix, celfix, spcfix, cylfix.

has_cdi_ja()

source code 

Returns True if CDi_ja is present. CDi_ja is an alternate specification of the linear transformation matrix, maintained for historical compatibility.

Matrix elements in the IRAF convention are equivalent to the product CDi_ja = CDELTia * PCi_ja, but the defaults differ from that of the PCi_ja matrix. If one or more CDi_ja keywords are present then all unspecified CDi_ja default to zero. If no CDi_ja (or CROTAia) keywords are present, then the header is assumed to be in PCi_ja form whether or not any PCi_ja keywords are present since this results in an interpretation of CDELTia consistent with the original FITS specification.

While CDi_ja may not formally co-exist with PCi_ja, it may co-exist with CDELTia and CROTAia which are to be ignored.

Returns: bool

See Also: cd for more information.

has_crotaia()

source code 

Returns True if CROTAia is present. CROTAia is an alternate specification of the linear transformation matrix, maintained for historical compatibility.

In the AIPS convention, CROTAia may only be associated with the latitude axis of a celestial axis pair. It specifies a rotation in the image plane that is applied AFTER the CDELTia; any other CROTAia keywords are ignored.

CROTAia may not formally co-exist with PCi_ja. CROTAia and CDELTia may formally co-exist with CDi_ja but if so are to be ignored.

Returns: bool

See Also: cd for more information.

has_pci_ja()

source code 

Returns True if PCi_ja is present. PCi_ja is the recommended way to specify the linear transformation matrix.

Returns: bool

See Also: cd for more information.

mix(mixpix, mixcel, vspan, vstep, viter, world, pixcrd)

source code 

Given either the celestial longitude or latitude plus an element of the pixel coordinate, solves for the remaining elements by iterating on the unknown celestial coordinate element using s2p.

Parameters:
  • mixpix (int) - Which element on the pixel coordinate is given.
  • mixcel (int) - Which element of the celestial coordinate is given. If 1, celestial longitude is given in world[self.lng], latitude returned in world[self.lat]. If 2, celestial latitude is given in world[self.lat], longitude returned in world[self.lng]
  • vspan (sequence of two floats) - Solution interval for the celestial coordinate, in degrees. The ordering of the two limits is irrelevant. Longitude ranges may be specified with any convenient normalization, for example (-120,+120) is the same as (240,480), except that the solution will be returned with the same normalization, i.e. lie within the interval specified.
  • vstep (float) - Step size for solution search, in degrees. If 0, a sensible, although perhaps non-optimal default will be used.
  • viter (int) - If a solution is not found then the step size will be halved and the search recommenced. viter controls how many times the step size is halved. The allowed range is 5 - 10.
  • world (array[naxis] of double) - World coordinate elements. world[self.lng] and world[self.lat] are the celestial longitude and latitude, in degrees. Which is given and which returned depends on the value of mixcel. All other elements are given. The results will be written to this array in-place.
  • pixcrd (array[naxis] of double) - Pixel coordinate. The element indicated by mixpix is given and the remaining elements will be written in-place. Pixel coordinates are zero-based.
Returns: dict
A dictionary with the following keys:
  • phi (array[naxis] of double)
  • theta (array[naxis] of double)
    • Longitude and latitude in the native coordinate system of the projection, in degrees.
  • imgcrd (array[naxis] of double)
    • Image coordinate elements. imgcrd[self.lng] and imgcrd[self.lat] are the projected x- and y-coordinates, in "degrees".
Raises:
  • MemoryError - Memory allocation failed.
  • SingularMatrixError - Linear transformation matrix is singular.
  • InconsistentAxisTypesError - Inconsistent or unrecognized coordinate axis types.
  • ValueError - Invalid parameter value.
  • InvalidTransformError - Invalid coordinate transformation parameters.
  • InvalidTransformError - Ill-conditioned coordinate transformation parameters.
  • InvalidCoordinateError - Invalid world coordinate.
  • NoSolutionError - No solution found in the specified interval.

p2s(pixcrd)

source code 

Converts pixel to world coordinates.

Parameters:
  • pixcrd (numpy array[ncoord][nelem] of double) - Array of pixel coordinates. Pixel coordinates are zero-based.
Returns: dict
A dictionary with the following keys:
  • imgcrd (array[ncoord][nelem] of double)
    • Array of intermediate world coordinates. For celestial axes, imgcrd[][self.lng] and imgcrd[][self.lat] are the projected x-, and y-coordinates, in decimal degrees. For spectral axes, imgcrd[][self.spec] is the intermediate spectral coordinate, in SI units.
  • phi (array[ncoord] of double)
  • theta (array[ncoord] of double)
    • Longitude and latitude in the native coordinate system of the projection, in degrees.
  • world (array[ncoord][nelem] of double)
    • Array of world coordinates. For celestial axes, world[][self.lng] and world[][self.lat] are the celestial longitude and latitude, in decimal degrees. For spectral axes, world[][self.spec] is the intermediate spectral coordinate, in SI units.
  • stat (array[ncoord] of int)
    • Status return value for each coordinate. 0 for success, 1 for invalid pixel coordinate.
Raises:
  • MemoryError - Memory allocation failed.
  • SingularMatrixError - Linear transformation matrix is singular.
  • InconsistentAxisTypesError - Inconsistent or unrecognized coordinate axis types.
  • ValueError - Invalid parameter value.
  • ValueError - x- and y-coordinate arrays are not the same size.
  • InvalidTransformError - Invalid coordinate transformation parameters.
  • InvalidTransformError - Ill-conditioned coordinate transformation parameters.

print_contents()

source code 

Print the contents of the WCS object to stdout. Probably only useful for debugging purposes, and may be removed in the future.

s2p(world)

source code 

Transforms world coordinates to pixel coordinates.

Parameters:
  • world (array[ncoord][nelem] of double) - Array of world coordinates, in decimal degrees.
Returns: dict
A dictionary with the following keys:
  • phi (array[ncoord] of double)
  • theta (array[ncoord] of double)
    • Longitude and latitude in the native coordinate system of the projection, in degrees.
  • imgcrd (array[ncoord][nelem] of double)
    • Array of intermediate world coordinates. For celestial axes, imgcrd[][self.lng] and imgcrd[][self.lat] are the projected x-, and y-coordinates, in "degrees". For quadcube projections with a CUBEFACE axis, the face number is also returned in imgcrd[][self.cubeface]. For spectral axes, imgcrd[][self.spec] is the intermediate spectral coordinate, in SI units.
  • pixcrd (array[ncoord][nelem] of double)
    • Array of pixel coordinates. Pixel coordinates are zero-based.
  • stat (array[ncoord] of int)
    • Status return value for each coordinate. 0 for success, 1 for invalid pixel coordinate.
Raises:
  • MemoryError - Memory allocation failed.
  • SingularMatrixError - Linear transformation matrix is singular.
  • InconsistentAxisTypesError - Inconsistent or unrecognized coordinate axis types.
  • ValueError - Invalid parameter value.
  • InvalidTransformError - Invalid coordinate transformation parameters.
  • InvalidTransformError - Ill-conditioned coordinate transformation parameters.

set()

source code 

Sets up a WCS object for use according to information supplied within it.

set recognizes the NCP projection and converts it to the equivalent SIN projection and it also recognizes GLS as a synonym for SFL. It does alias translation for the AIPS spectral types (FREQ-LSR, FELO-HEL, etc.) but without changing the input header keywords.

Note that this routine need not be called directly; it will be invoked by p2s and s2p if necessary.

set strips off trailing blanks in all string members.

Raises:
  • MemoryError - Memory allocation failed.
  • SingularMatrixError - Linear transformation matrix is singular.
  • InconsistentAxisTypesError - Inconsistent or unrecognized coordinate axis types.
  • ValueError - Invalid parameter value.
  • InvalidTransformError - Invalid coordinate transformation parameters.
  • InvalidTransformError - Ill-conditioned coordinate transformation parameters.

spcfix()

source code 

Translates AIPS-convention spectral coordinate types. {FREQ,VELO,FELO}-{OBS,HEL,LSR} (e.g. FREQ-LSR, VELO-OBS, FELO-HEL)

Returns: int
0 for success; -1 if no change required.

sptr(ctype, i=-1)

source code 

Translates the spectral axis in a WCS object. For example, a FREQ axis may be translated into ZOPT-F2W and vice versa.

Parameters:
  • ctype (string) - Required spectral CTYPEia. Wildcarding may be used, i.e. if the final three characters are specified as "???", or if just the eighth character is specified as "?", the correct algorithm code will be substituted and returned.
  • i (int) - Index of the spectral axis (0-rel). If i < 0 (or not provided), it will be set to the first spectral axis identified from the CTYPE keyvalues in the FITS header.
Raises:
  • MemoryError - Memory allocation failed.
  • SingularMatrixError - Linear transformation matrix is singular.
  • InconsistentAxisTypesError - Inconsistent or unrecognized coordinate axis types.
  • ValueError - Invalid parameter value.
  • InvalidTransformError - Invalid coordinate transformation parameters.
  • InvalidTransformError - Ill-conditioned coordinate transformation parameters.
  • InvalidSubimageSpecificationError - Invalid subimage specification (no spectral axis).

unitfix(translate_units='')

source code 

Translates non-standard CUNITia keyvalues. For example, DEG -> deg, also stripping off unnecessary whitespace.

Parameters:
  • translate_units - Do potentially unsafe translations of non-standard unit strings.

    Although "S" is commonly used to represent seconds, its recognizes "S" formally as Siemens, however rarely that may be translation to "s" is potentially unsafe since the standard used. The same applies to "H" for hours (Henry), and "D" for days (Debye).

    This string controls what to do in such cases, and is case-insensitive.

    • If the string contains "s", translate "S" to "s".
    • If the string contains "h", translate "H" to "h".
    • If the string contains "d", translate "D" to "d".

    Thus '' doesn't do any unsafe translations, whereas 'shd' does all of them.

Returns: int
0 for success; -1 if no change required.

Property Details

alt

Character code for alternate coordinate descriptions. For example, the "a" in keyword names such as CTYPEia). This is blank for the primary coordinate description, or one of the 26 upper-case letters, A-Z.

Type:
string

cd

CDi_ja linear transformation matrix.

For historical compatibility, two alternate specifications of the CDi_ja and CROTAia keywords. Although these may not formally co-exist with PCi_ja, the approach here is simply to ignore them if given in conjunction with PCi_ja.

has_pci_ja, has_cdi_ja and has_crotaia can be used to determine which of these alternatives are present in the header.

CDi_ja and CROTAia keywords, if found, are to be stored in the cd and crota arrays which are dimensioned similarly to pc and cdelt.

These alternate specifications of the linear transformation matrix are translated immediately to PCi_ja by set and are nowhere visible to the lower-level routines. In particular, set resets cdelt to unity if CDi_ja is present (and no PCi_ja). If no CROTAia is associated with the latitude axis, set reverts to a unity PCi_ja matrix.

Type:
array[2][2] of double

cdelt

Coordinate increments (CDELTia) for each coord axis.

Any undefined elements are returned as NaN.

Type:
array[naxis] of double

colnum

Where the coordinate representation is associated with an image-array column in a FITS binary table, this variable may be used to record the relevant column number.

It should be set to zero for an image header or pixel list.

crder

The random error in each coordinate axis, CRDERia.

Any undefined elements are returned as NaN.

Type:
array[naxis] of double

See Also: csyer

crota

CROTAia keyvalues for each coordinate axis.

CROTAia is an alternate specification of the linear transformation matrix, maintained for historical compatibility.

Type:
array[2][2] of double

See Also: cd for more information.

crpix

Coordinate reference pixels (CRPIXja) for each pixel axis.

Pixel coordinates are zero-based.

Type:
array[naxis] of double

csyer

The systematic error in the coordinate value axes, CSYERia.

Any undefined elements are returned as NaN.

Type:
array[naxis] of double

See Also: crder

ctype

List of CTYPEia keyvalues.

The ctype keyword values must be in upper case and there must be zero or one pair of matched celestial axis types, and zero or one spectral axis.

Type:
list of strings

cubeface

Index into the pixcrd (pixel coordinate) array for the CUBEFACE axis. This is used for quadcube projections where the cube faces are stored on a separate axis.

The quadcube projections (TSC, CSC, QSC) may be represented in FITS in either of two ways:

  • The six faces may be laid out in one plane and numbered as follows:
                                    0
    
                           4  3  2  1  4  3  2
    
                                    5
    

    Faces 2, 3 and 4 may appear on one side or the other (or both). The world-to-pixel routines map faces 2, 3 and 4 to the left but the pixel-to-world routines accept them on either side.

  • The "COBE" convention in which the six faces are stored in a three-dimensional structure using a "CUBEFACE" axis indexed from 0 to 5 as above.

These routines support both methods; set determines which is being used by the presence or absence of a CUBEFACE axis in ctype. p2s and s2p translate the CUBEFACE axis representation to the single plane representation understood by the lower-level projection routines.

Type:
int

cunit

List of CUNITia keyvalues which define the units of measurement of the CRVALia, CDELTia and CDi_ja keywords.

As CUNITia is an optional header keyword, cunit may be left blank but otherwise is expected to contain a standard units specification as defined by WCS Paper I. Utility function wcsutrn(), (not currently wrapped for Python) is available to translate commonly used non-standard units specifications but this must be done as a separate step before invoking set.

For celestial axes, if cunit is not blank, set uses wcsunits to parse it and scale cdelt, crval, and cd to decimal degrees. It then resets cunit to "deg".

For spectral axes, if cunit is not blank, set uses wcsunits to parse it and scale cdelt, crval, and cd to SI units. It then resets cunit accordingly.

set ignores cunit for other coordinate types; cunit may be used to label coordinate values.

Type:
list of strings

dateavg

Representative mid-point of the date of observation in ISO format, yyyy-mm-ddThh:mm:ss.

Type:
string

See Also: dateobs

dateobs

Start of the date of observation in ISO format, yyyy-mm-ddThh:mm:ss.

Type:
string

See Also: dateavg

equinox

The equinox associated with dynamical equatorial or ecliptic coordinate systems, EQUINOXa (or EPOCH in older headers). Not applicable to ICRS equatorial or ecliptic coordinates.

If the equinox is undefined, NaN is returned.

Type:
float

lat

The index into the world coordinate array containing latitude values.

Type:
int

See Also: lng

latpole

The native latitude of the celestial pole, LATPOLEa (deg).

Type:
float

See Also: lonpole

lng

The index into the world coordinate array containing longitude values.

Type:
int

See Also: lat

lonpole

The native longitude of the celestial pole, LONPOLEa (deg).

Type:
float

See Also: latpole

mjdavg

Modified Julian Date (MJD = JD - 2400000.5), MJD-AVG, corresponding to DATE-AVG.

If mjdavg is undefined, NaN is returned.

Type:
float

See Also: mjdobs

mjdobs

Modified Julian Date (MJD = JD - 2400000.5), MJD-OBS, corresponding to DATE-OBS.

If mjdobs is undefined, NaN is returned.

Type:
float

See Also: mjdavg

obsgeo

Location of the observer in a standard terrestrial reference frame, OBSGEO-X, OBSGEO-Y, OBSGEO-Z (in metres).

Any undefined elements are returned as NaN.

Type:
array[3] of double

pc

The PCi_ja (pixel coordinate) transformation matrix. The order is:

 [[PC1_1, PC1_2],
  [PC2_1, PC2_2]]
Type:
array[2][2] of double

ps

PSi_ma keywords for each i and m. Returned as a list of tuples of the form (i, m, value):

  • i: axis number, as in PSi_ma, (i.e. 1-relative)
  • m: parameter number, as in PSi_ma, (i.e. 0-relative)
  • value: parameter value (as a string)
Type:
list of tuples

pv

PVi_ma keywords for each i and m. Returned as a list of tuples of the form (i, m, value):

  • i: axis number, as in PVi_ma, (i.e. 1-relative)
  • m: parameter number, as in PVi_ma, (i.e. 0-relative)
  • value: parameter value (as a string)
Type:
list of tuples

restfrq

Rest frequency (Hz) from RESTFRQa.

If restfrq is undefined, NaN is returned.

Type:
float

See Also: restwav

restwav

Rest wavelength (m) from RESTWAVa.

If restwav is undefined, NaN is returned.

Type:
float

See Also: restfrq

specsys

Spectral reference frame (standard of rest), SPECSYSa

Type:
string

See Also: ssysobs, velosys.

ssysobs

The actual spectral reference frame in which there is no differential variation in the spectral coordinate across the field-of-view, SSYSOBSa.

Type:
string

See Also: specsys, velosys

ssyssrc

The spectral reference frame (standard of rest) in which the redshift was measured, SSYSSRCa.

Type:
string

See Also: zsource

velangl

The angle in degrees that should be used to decompose an observed velocity into radial and transverse components.

If velangl is undefined, NaN is returned.

Type:
float

velosys

The relative radial velocity (m/s) between the observer and the selected standard of rest in the direction of the celestial reference coordinate, VELOSYSa.

If velosys is undefined, NaN is returned.

Type:
float

See Also: specsys, ssysobs

zsource

The redshift, ZSOURCEa, of the source.

If zsource is undefined, NaN is returned.

Type:
float

See Also: ssyssrc