Wcsprm#

class astropy.wcs.Wcsprm(header=None, key=' ', relax=False, naxis=2, keysel=0, colsel=None)#

Bases: object

Wcsprm performs the core WCS transformations.

Note

The members of this object correspond roughly to the key/value pairs in the FITS header. However, they are adjusted and normalized in a number of ways that make performing the WCS transformation easier. Therefore, they can not be relied upon to get the original values in the header. For that, use astropy.io.fits.Header directly.

The FITS header parsing enforces correct FITS “keyword = value” syntax with regard to the equals sign occurring in columns 9 and 10. However, it does recognize free-format character (NOST 100-2.0, Sect. 5.2.1), integer (Sect. 5.2.3), and floating-point values (Sect. 5.2.4) for all keywords.

Warning

Many of the attributes of this class require additional processing when modifying underlying C structure. When needed, this additional processing is implemented in attribute setters. Therefore, for mutable attributes, one should always set the attribute rather than a slice of its current value (or its individual elements) since the latter may lead the class instance to be in an invalid state. For example, attribute crpix of a 2D WCS’ Wcsprm object wcs should be set as wcs.crpix = [crpix1, crpix2] instead of wcs.crpix[0] = crpix1; wcs.crpix[1] = crpix2].

Parameters:
headerHeader, str, or None.

If None, the object will be initialized to default values.

keystr, optional

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". (key may only be provided if header is also provided.)

relaxbool or int, optional

Degree of permissiveness:

  • False: Recognize only FITS keywords defined by the published WCS standard.

  • True: Admit all recognized informal extensions of the WCS standard.

  • int: a bit field selecting specific extensions to accept. See Header-reading relaxation constants for details.

naxisint, optional

The number of world coordinates axes for the object. (naxis may only be provided if header is None.)

keyselsequence of flag bits, optional

Vector of flag bits that may be used to restrict the keyword types considered:

  • WCSHDR_IMGHEAD: Image header keywords.

  • WCSHDR_BIMGARR: Binary table image array.

  • WCSHDR_PIXLIST: Pixel list keywords.

If zero, there is no restriction. If -1, the underlying wcslib function wcspih() is called, rather than wcstbh().

colselsequence of int

A sequence of table column numbers used to restrict the keywords considered. None indicates no restriction.

Raises:
MemoryError

Memory allocation failed.

ValueError

Invalid key.

KeyError

Key not found in FITS header.

Attributes Summary

alt

str Character code for alternate coordinate descriptions.

aux

Auxprm Auxiliary coordinate system information of a specialist nature.

axis_types

int array[naxis] An array of four-digit type codes for each axis.

bepoch

double Equivalent to DATE-OBS.

cd

double array[naxis][naxis] The CDi_ja linear transformation matrix.

cdelt

double array[naxis] Coordinate increments (CDELTia) for each coord axis.

cel

Celprm Information required to transform celestial coordinates.

cel_offset

boolean Is there an offset?

cname

list of strings A list of the coordinate axis names, from CNAMEia.

colax

int array[naxis] An array recording the column numbers for each axis in a pixel list.

colnum

int Column of FITS binary table associated with this WCS.

cperi

double array[naxis] period of a phase axis, CPERIia.

crder

double array[naxis] The random error in each coordinate axis, CRDERia.

crota

double array[naxis] CROTAia keyvalues for each coordinate axis.

crpix

double array[naxis] Coordinate reference pixels (CRPIXja) for each pixel axis.

crval

double array[naxis] Coordinate reference values (CRVALia) for each coordinate axis.

csyer

double array[naxis] The systematic error in the coordinate value axes, CSYERia.

ctype

list of strings[naxis] List of CTYPEia keyvalues.

cubeface

int Index into the pixcrd (pixel coordinate) array for the CUBEFACE axis.

cunit

list of astropy.UnitBase[naxis] List of CUNITia keyvalues as astropy.units.UnitBase instances.

czphs

double array[naxis] The time at the zero point of a phase axis, CSPHSia.

dateavg

string Representative mid-point of the date of observation.

datebeg

string Date at the start of the observation.

dateend

string Date at the end of the observation.

dateobs

string Start of the date of observation.

dateref

string Date of a reference epoch relative to which other time measurements refer.

equinox

double The equinox associated with dynamical equatorial or ecliptic coordinate systems.

imgpix_matrix

double array[2][2] (read-only) Inverse of the CDELT or PC matrix.

jepoch

double Equivalent to DATE-OBS.

lat

int (read-only) The index into the world coord array containing latitude values.

latpole

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

lattyp

string (read-only) Celestial axis type for latitude.

lng

int (read-only) The index into the world coord array containing longitude values.

lngtyp

string (read-only) Celestial axis type for longitude.

lonpole

double The native longitude of the celestial pole.

mjdavg

double Modified Julian Date corresponding to DATE-AVG.

mjdbeg

double Modified Julian Date corresponding to DATE-BEG.

mjdend

double Modified Julian Date corresponding to DATE-END.

mjdobs

double Modified Julian Date corresponding to DATE-OBS.

mjdref

double Modified Julian Date corresponding to DATE-REF.

name

string The name given to the coordinate representation WCSNAMEa.

naxis

int (read-only) The number of axes (pixel and coordinate).

obsgeo

double array[3] Location of the observer in a standard terrestrial reference frame.

obsorbit

string URI, URL, or name of an orbit ephemeris file giving spacecraft coordinates relating to TREFPOS.

pc

double array[naxis][naxis] The PCi_ja (pixel coordinate) transformation matrix.

phi0

double The native latitude of the fiducial point.

piximg_matrix

double array[2][2] (read-only) Matrix containing the product of the CDELTia diagonal matrix and the PCi_ja matrix.

plephem

string The Solar System ephemeris used for calculating a pathlength delay.

radesys

string The equatorial or ecliptic coordinate system type, RADESYSa.

restfrq

double Rest frequency (Hz) from RESTFRQa.

restwav

double Rest wavelength (m) from RESTWAVa.

spec

int (read-only) The index containing the spectral axis values.

specsys

string Spectral reference frame (standard of rest), SPECSYSa.

ssysobs

string Spectral reference frame.

ssyssrc

string Spectral reference frame for redshift.

tab

list of Tabprm Tabular coordinate objects.

telapse

double equivalent to the elapsed time between DATE-BEG and DATE-END, in units of TIMEUNIT.

theta0

double The native longitude of the fiducial point.

timedel

double the resolution of the time stamps.

timeoffs

double Time offset, which may be used, for example, to provide a uniform clock correction

timepixr

double relative position of the time stamps in binned time intervals, a value between 0.0 and 1.0.

timesys

string Time scale (UTC, TAI, etc.) in which all other time-related auxiliary header values are recorded.

timeunit

string Time units in which the following header values are expressed: TSTART, TSTOP, TIMEOFFS, TIMSYER, TIMRDER, TIMEDEL.

timrder

double the accuracy of time stamps relative to each other, in units of TIMEUNIT.

timsyer

double the absolute error of the time values, in units of TIMEUNIT.

trefdir

string Reference direction used in calculating a pathlength delay.

trefpos

string Location in space where the recorded time is valid.

tstart

double equivalent to DATE-BEG expressed as a time in units of TIMEUNIT relative to DATEREF+TIMEOFFS.

tstop

double equivalent to DATE-END expressed as a time in units of TIMEUNIT relative to DATEREF+TIMEOFFS.

velangl

double Velocity angle.

velosys

double Relative radial velocity.

velref

int AIPS velocity code.

wtb

list of Wtbarr objects to construct coordinate lookup tables from BINTABLE.

xposure

double effective exposure time in units of TIMEUNIT.

zsource

double The redshift, ZSOURCEa, of the source.

Methods Summary

bounds_check(pix2world, world2pix)

Enable/disable bounds checking.

cdfix()

Fix erroneously omitted CDi_ja keywords.

celfix

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

compare(other[, cmp, tolerance])

Compare two Wcsprm objects for equality.

cylfix()

Fixes WCS keyvalues for malformed cylindrical projections.

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.

fix([translate_units, naxis])

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

get_cdelt()

Coordinate increments (CDELTia) for each coord axis as double array[naxis].

get_pc()

Returns the PC matrix in read-only form as double array[naxis][naxis].

get_ps()

Returns PSi_ma keywords for each i and m as list of tuples.

get_pv()

Returns PVi_ma keywords for each i and m as list of tuples.

has_cd()

Returns True if CDi_ja is present.

has_cdi_ja()

Alias for has_cd.

has_crota()

Returns True if CROTAia is present.

has_crotaia()

Alias for has_crota.

has_pc()

Returns True if PCi_ja is present.

has_pci_ja()

Alias for has_pc.

is_unity()

Returns True if the linear transformation matrix (cd) is unity.

mix(mixpix, mixcel, vspan, vstep, viter, ...)

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.

p2s(pixcrd, origin)

Converts pixel to world coordinates.

print_contents()

Print the contents of the Wcsprm object to stdout.

s2p(world, origin)

Transforms world coordinates to pixel coordinates.

set()

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

set_ps(ps)

Sets PSi_ma keywords for each i and m.

set_pv(pv)

Sets PVi_ma keywords for each i and m.

spcfix()

Translates AIPS-convention spectral coordinate types.

sptr(ctype[, i])

Translates the spectral axis in a WCS object.

sub(axes)

Extracts the coordinate description for a subimage from a WCS object.

to_header([relax])

to_header translates a WCS object into a FITS header.

unitfix([translate_units])

Translates non-standard CUNITia keyvalues.

Attributes Documentation

alt#

str Character code for alternate coordinate descriptions.

For example, the "a" in keyword names such as CTYPEia. This is a space character for the primary coordinate description, or one of the 26 upper-case letters, A-Z.

aux#

Auxprm Auxiliary coordinate system information of a specialist nature.

axis_types#

int array[naxis] An array of four-digit type codes for each axis.

  • First digit (i.e. 1000s):

    • 0: Non-specific coordinate type.

    • 1: Stokes coordinate.

    • 2: Celestial coordinate (including CUBEFACE).

    • 3: Spectral coordinate.

  • Second digit (i.e. 100s):

    • 0: Linear axis.

    • 1: Quantized axis (STOKES, CUBEFACE).

    • 2: Non-linear celestial axis.

    • 3: Non-linear spectral axis.

    • 4: Logarithmic axis.

    • 5: Tabular axis.

  • Third digit (i.e. 10s):

    • 0: Group number, e.g. lookup table number

  • The fourth digit is used as a qualifier depending on the axis type.

    • For celestial axes:

      • 0: Longitude coordinate.

      • 1: Latitude coordinate.

      • 2: CUBEFACE number.

    • For lookup tables: the axis number in a multidimensional table.

CTYPEia in "4-3" form with unrecognized algorithm code will have its type set to -1 and generate an error.

bepoch#

double Equivalent to DATE-OBS.

Expressed as a Besselian epoch.

cd#

double array[naxis][naxis] The CDi_ja linear transformation matrix.

For historical compatibility, three alternate specifications of the linear transformations are available in wcslib. The canonical PCi_ja with CDELTia, CDi_ja, and the deprecated CROTAia keywords. Although the latter may not formally co-exist with PCi_ja, the approach here is simply to ignore them if given in conjunction with PCi_ja.

has_pc, has_cd and has_crota can be used to determine which of these alternatives are present in the header.

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.

cdelt#

double array[naxis] Coordinate increments (CDELTia) for each coord axis.

If a CDi_ja linear transformation matrix is present, a warning is raised and cdelt is ignored. The CDi_ja matrix may be deleted by:

del wcs.wcs.cd

An undefined value is represented by NaN.

cel#

Celprm Information required to transform celestial coordinates.

cel_offset#

boolean Is there an offset?

If True, an offset will be applied to (x, y) to force (x, y) = (0, 0) at the fiducial point, (phi_0, theta_0). Default is False.

cname#

list of strings A list of the coordinate axis names, from CNAMEia.

colax#

int array[naxis] An array recording the column numbers for each axis in a pixel list.

colnum#

int Column of FITS binary table associated with this WCS.

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

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

cperi#

double array[naxis] period of a phase axis, CPERIia.

An undefined value is represented by NaN.

crder#

double array[naxis] The random error in each coordinate axis, CRDERia.

An undefined value is represented by NaN.

crota#

double array[naxis] CROTAia keyvalues for each coordinate axis.

For historical compatibility, three alternate specifications of the linear transformations are available in wcslib. The canonical PCi_ja with CDELTia, CDi_ja, and the deprecated CROTAia keywords. Although the latter may not formally co-exist with PCi_ja, the approach here is simply to ignore them if given in conjunction with PCi_ja.

has_pc, has_cd and has_crota can be used to determine which of these alternatives are present in the header.

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.

crpix#

double array[naxis] Coordinate reference pixels (CRPIXja) for each pixel axis.

crval#

double array[naxis] Coordinate reference values (CRVALia) for each coordinate axis.

csyer#

double array[naxis] The systematic error in the coordinate value axes, CSYERia.

An undefined value is represented by NaN.

ctype#

list of strings[naxis] 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.

cubeface#

int 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.

cunit#

list of astropy.UnitBase[naxis] List of CUNITia keyvalues as astropy.units.UnitBase instances.

These 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. unitfix 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.

czphs#

double array[naxis] The time at the zero point of a phase axis, CSPHSia.

An undefined value is represented by NaN.

dateavg#

string Representative mid-point of the date of observation.

In ISO format, yyyy-mm-ddThh:mm:ss.

datebeg#

string Date at the start of the observation.

In ISO format, yyyy-mm-ddThh:mm:ss.

dateend#

string Date at the end of the observation.

In ISO format, yyyy-mm-ddThh:mm:ss.

dateobs#

string Start of the date of observation.

In ISO format, yyyy-mm-ddThh:mm:ss.

dateref#

string Date of a reference epoch relative to which other time measurements refer.

equinox#

double The equinox associated with dynamical equatorial or ecliptic coordinate systems.

EQUINOXa (or EPOCH in older headers). Not applicable to ICRS equatorial or ecliptic coordinates.

An undefined value is represented by NaN.

imgpix_matrix#

double array[2][2] (read-only) Inverse of the CDELT or PC matrix.

Inverse containing the product of the CDELTia diagonal matrix and the PCi_ja matrix.

jepoch#

double Equivalent to DATE-OBS.

Expressed as a Julian epoch.

lat#

int (read-only) The index into the world coord array containing latitude values.

latpole#

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

lattyp#

string (read-only) Celestial axis type for latitude.

For example, “RA”, “DEC”, “GLON”, “GLAT”, etc. extracted from “RA–“, “DEC-”, “GLON”, “GLAT”, etc. in the first four characters of CTYPEia but with trailing dashes removed.

lng#

int (read-only) The index into the world coord array containing longitude values.

lngtyp#

string (read-only) Celestial axis type for longitude.

For example, “RA”, “DEC”, “GLON”, “GLAT”, etc. extracted from “RA–“, “DEC-”, “GLON”, “GLAT”, etc. in the first four characters of CTYPEia but with trailing dashes removed.

lonpole#

double The native longitude of the celestial pole.

LONPOLEa (deg).

mjdavg#

double Modified Julian Date corresponding to DATE-AVG.

(MJD = JD - 2400000.5).

An undefined value is represented by NaN.

mjdbeg#

double Modified Julian Date corresponding to DATE-BEG.

(MJD = JD - 2400000.5).

An undefined value is represented by NaN.

mjdend#

double Modified Julian Date corresponding to DATE-END.

(MJD = JD - 2400000.5).

An undefined value is represented by NaN.

mjdobs#

double Modified Julian Date corresponding to DATE-OBS.

(MJD = JD - 2400000.5).

An undefined value is represented by NaN.

mjdref#

double Modified Julian Date corresponding to DATE-REF.

(MJD = JD - 2400000.5).

An undefined value is represented by NaN.

name#

string The name given to the coordinate representation WCSNAMEa.

naxis#

int (read-only) The number of axes (pixel and coordinate).

Given by the NAXIS or WCSAXESa keyvalues.

The number of coordinate axes is determined at parsing time, and can not be subsequently changed.

It is determined from the highest of the following:

  1. NAXIS

  2. WCSAXESa

  3. The highest axis number in any parameterized WCS keyword. The keyvalue, as well as the keyword, must be syntactically valid otherwise it will not be considered.

If none of these keyword types is present, i.e. if the header only contains auxiliary WCS keywords for a particular coordinate representation, then no coordinate description is constructed for it.

This value may differ for different coordinate representations of the same image.

obsgeo#

double array[3] Location of the observer in a standard terrestrial reference frame.

OBSGEO-X, OBSGEO-Y, OBSGEO-Z (in meters).

An undefined value is represented by NaN.

obsorbit#

string URI, URL, or name of an orbit ephemeris file giving spacecraft coordinates relating to TREFPOS.

pc#

double array[naxis][naxis] The PCi_ja (pixel coordinate) transformation matrix.

The order is:

[[PC1_1, PC1_2],
 [PC2_1, PC2_2]]

For historical compatibility, three alternate specifications of the linear transformations are available in wcslib. The canonical PCi_ja with CDELTia, CDi_ja, and the deprecated CROTAia keywords. Although the latter may not formally co-exist with PCi_ja, the approach here is simply to ignore them if given in conjunction with PCi_ja.

has_pc, has_cd and has_crota can be used to determine which of these alternatives are present in the header.

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.

phi0#

double The native latitude of the fiducial point.

The point whose celestial coordinates are given in ref[1:2]. If undefined (NaN) the initialization routine, set, will set this to a projection-specific default.

piximg_matrix#

double array[2][2] (read-only) Matrix containing the product of the CDELTia diagonal matrix and the PCi_ja matrix.

plephem#

string The Solar System ephemeris used for calculating a pathlength delay.

radesys#

string The equatorial or ecliptic coordinate system type, RADESYSa.

restfrq#

double Rest frequency (Hz) from RESTFRQa.

An undefined value is represented by NaN.

restwav#

double Rest wavelength (m) from RESTWAVa.

An undefined value is represented by NaN.

spec#

int (read-only) The index containing the spectral axis values.

specsys#

string Spectral reference frame (standard of rest), SPECSYSa.

ssysobs#

string Spectral reference frame.

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

ssyssrc#

string Spectral reference frame for redshift.

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

tab#

list of Tabprm Tabular coordinate objects.

A list of tabular coordinate objects associated with this WCS.

telapse#

double equivalent to the elapsed time between DATE-BEG and DATE-END, in units of TIMEUNIT.

theta0#

double The native longitude of the fiducial point.

The point whose celestial coordinates are given in ref[1:2]. If undefined (NaN) the initialization routine, set, will set this to a projection-specific default.

timedel#

double the resolution of the time stamps.

timeoffs#
double Time offset, which may be used, for example, to provide a uniform clock correction

for times referenced to DATEREF.

timepixr#

double relative position of the time stamps in binned time intervals, a value between 0.0 and 1.0.

timesys#

string Time scale (UTC, TAI, etc.) in which all other time-related auxiliary header values are recorded. Also defines the time scale for an image axis with CTYPEia set to ‘TIME’.

timeunit#

string Time units in which the following header values are expressed: TSTART, TSTOP, TIMEOFFS, TIMSYER, TIMRDER, TIMEDEL.

It also provides the default value for CUNITia for time axes.

timrder#

double the accuracy of time stamps relative to each other, in units of TIMEUNIT.

timsyer#

double the absolute error of the time values, in units of TIMEUNIT.

trefdir#

string Reference direction used in calculating a pathlength delay.

trefpos#

string Location in space where the recorded time is valid.

tstart#

double equivalent to DATE-BEG expressed as a time in units of TIMEUNIT relative to DATEREF+TIMEOFFS.

tstop#

double equivalent to DATE-END expressed as a time in units of TIMEUNIT relative to DATEREF+TIMEOFFS.

velangl#

double Velocity angle.

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

An undefined value is represented by NaN.

velosys#

double Relative radial velocity.

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

An undefined value is represented by NaN.

velref#

int AIPS velocity code.

From VELREF keyword.

wtb#

list of Wtbarr objects to construct coordinate lookup tables from BINTABLE.

xposure#

double effective exposure time in units of TIMEUNIT.

zsource#

double The redshift, ZSOURCEa, of the source.

An undefined value is represented by NaN.

Methods Documentation

bounds_check(pix2world, world2pix)#

Enable/disable bounds checking.

Parameters:
pix2worldbool, optional

When True, enable bounds checking for the pixel-to-world (p2x) transformations. Default is True.

world2pixbool, optional

When True, enable bounds checking for the world-to-pixel (s2x) transformations. Default is True.

Notes

Note that by default (without calling bounds_check) strict bounds checking is enabled.

cdfix()#

Fix erroneously omitted CDi_ja keywords.

Sets the diagonal element of the CDi_ja matrix to unity if all CDi_ja keywords associated with a given axis were omitted. According to Paper I, if any CDi_ja keywords at all are given in a FITS header then those not given default to zero. This results in a singular matrix with an intersecting row and column of zeros.

Returns:
successint

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

celfix()#

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

Returns:
successint

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

compare(other, cmp=0, tolerance=0.0)#

Compare two Wcsprm objects for equality.

Parameters:
otherWcsprm

The other Wcsprm object to compare to.

cmpint, optional

A bit field controlling the strictness of the comparison. When 0, (the default), all fields must be identical.

The following constants, defined in the astropy.wcs module, may be or’ed together to loosen the comparison.

  • WCSCOMPARE_ANCILLARY: Ignores ancillary keywords that don’t change the WCS transformation, such as XPOSURE or EQUINOX. Note that this also ignores DATE-OBS, which does change the WCS transformation in some cases.

  • WCSCOMPARE_TILING: Ignore integral differences in CRPIXja. This is the ‘tiling’ condition, where two WCSes cover different regions of the same map projection and align on the same map grid.

  • WCSCOMPARE_CRPIX: Ignore any differences at all in CRPIXja. The two WCSes cover different regions of the same map projection but may not align on the same grid map. Overrides WCSCOMPARE_TILING.

tolerancefloat, optional

The amount of tolerance required. For example, for a value of 1e-6, all floating-point values in the objects must be equal to the first 6 decimal places. The default value of 0.0 implies exact equality.

Returns:
equalbool
cylfix()#

Fixes WCS keyvalues for malformed cylindrical projections.

Returns:
successint

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

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.

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:
successint

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

fix(translate_units='', naxis=0)#

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

Parameters:
translate_unitsstr, optional

Specify which potentially unsafe translations of non-standard unit strings to perform. By default, performs all.

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.

naxisint array, optional

Image axis lengths. If this array is set to zero or None, then cylfix will not be invoked.

Returns:
statusdict

Returns a dictionary containing the following keys, each referring to a status string for each of the sub-fix functions that were called:

get_cdelt() numpy.ndarray#

Coordinate increments (CDELTia) for each coord axis as double array[naxis].

Returns the CDELT offsets in read-only form. Unlike the cdelt property, this works even when the header specifies the linear transformation matrix in one of the alternative CDi_ja or CROTAia forms. This is useful when you want access to the linear transformation matrix, but don’t care how it was specified in the header.

get_pc() numpy.ndarray#

Returns the PC matrix in read-only form as double array[naxis][naxis]. Unlike the pc property, this works even when the header specifies the linear transformation matrix in one of the alternative CDi_ja or CROTAia forms. This is useful when you want access to the linear transformation matrix, but don’t care how it was specified in the header.

get_ps() list#

Returns PSi_ma keywords for each i and m as list of tuples.

Returns:
pslist

Returned as a list of tuples of the form (i, m, value):

  • i: int. Axis number, as in PSi_ma, (i.e. 1-relative)

  • m: int. Parameter number, as in PSi_ma, (i.e. 0-relative)

  • value: string. Parameter value.

See also

astropy.wcs.Wcsprm.set_ps

Set PSi_ma values

get_pv() list#

Returns PVi_ma keywords for each i and m as list of tuples.

Returns:
sequence of tuple

Returned as a list of tuples of the form (i, m, value):

  • i: int. Axis number, as in PVi_ma, (i.e. 1-relative)

  • m: int. Parameter number, as in PVi_ma, (i.e. 0-relative)

  • value: string. Parameter value.

See also

astropy.wcs.Wcsprm.set_pv

Set PVi_ma values

Notes

Note that, if they were not given, set resets the entries for PVi_1a, PVi_2a, PVi_3a, and PVi_4a for longitude axis i to match (phi_0, theta_0), the native longitude and latitude of the reference point given by LONPOLEa and LATPOLEa.

has_cd() bool#

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.

See also

astropy.wcs.Wcsprm.cd

Get the raw CDi_ja values.

has_cdi_ja() bool#

Alias for has_cd. Maintained for backward compatibility.

has_crota() bool#

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.

See also

astropy.wcs.Wcsprm.crota

Get the raw CROTAia values

has_crotaia() bool#

Alias for has_crota. Maintained for backward compatibility.

has_pc() bool#

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

See also

astropy.wcs.Wcsprm.pc

Get the raw PCi_ja values

has_pci_ja() bool#

Alias for has_pc. Maintained for backward compatibility.

is_unity() bool#

Returns True if the linear transformation matrix (cd) is unity.

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

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:
mixpixint

Which element on the pixel coordinate is given.

mixcelint

Which element of the celestial coordinate is given. If mixcel = 1, celestial longitude is given in world[self.lng], latitude returned in world[self.lat]. If mixcel = 2, celestial latitude is given in world[self.lat], longitude returned in world[self.lng].

vspan(float, float)

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.

vstepfloat

Step size for solution search, in degrees. If 0, a sensible, although perhaps non-optimal default will be used.

viterint

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.

worldndarray

World coordinate elements as double array[naxis]. 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.

pixcrdndarray

Pixel coordinates as double array[naxis]. The element indicated by mixpix is given and the remaining elements will be written in-place.

originint

Specifies the origin of pixel values. The Fortran and FITS standards use an origin of 1. Numpy and C use array indexing with origin at 0.

Returns:
resultdict

Returns a dictionary with the following keys:

  • phi (double array[naxis])

  • theta (double array[naxis])

    • Longitude and latitude in the native coordinate system of the projection, in degrees.

  • imgcrd (double array[naxis])

    • Image coordinate elements. imgcrd[self.lng] and imgcrd[self.lat] are the projected x- and y-coordinates, in decimal degrees.

  • world (double array[naxis])

    • Another reference to the world argument passed in.

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.

See also

astropy.wcs.Wcsprm.lat, astropy.wcs.Wcsprm.lng

Get the axes numbers for latitude and longitude

Notes

Initially, the specified solution interval is checked to see if it’s a “crossing” interval. If it isn’t, a search is made for a crossing solution by iterating on the unknown celestial coordinate starting at the upper limit of the solution interval and decrementing by the specified step size. A crossing is indicated if the trial value of the pixel coordinate steps through the value specified. If a crossing interval is found then the solution is determined by a modified form of “regula falsi” division of the crossing interval. If no crossing interval was found within the specified solution interval then a search is made for a “non-crossing” solution as may arise from a point of tangency. The process is complicated by having to make allowance for the discontinuities that occur in all map projections.

Once one solution has been determined others may be found by subsequent invocations of mix with suitably restricted solution intervals.

Note the circumstance that arises when the solution point lies at a native pole of a projection in which the pole is represented as a finite curve, for example the zenithals and conics. In such cases two or more valid solutions may exist but mix only ever returns one.

Because of its generality, mix is very compute-intensive. For compute-limited applications, more efficient special-case solvers could be written for simple projections, for example non-oblique cylindrical projections.

p2s(pixcrd, origin)#

Converts pixel to world coordinates.

Parameters:
pixcrdndarray

Array of pixel coordinates as double array[ncoord][nelem].

originint

Specifies the origin of pixel values. The Fortran and FITS standards use an origin of 1. Numpy and C use array indexing with origin at 0.

Returns:
resultdict

Returns a dictionary with the following keys:

  • imgcrd: ndarray

    • Array of intermediate world coordinates as double array[ncoord][nelem]. For celestial axes, imgcrd[][self.lng] and imgcrd[][self.lat] are the projected x-, and y-coordinates, in pseudo degrees. For spectral axes, imgcrd[][self.spec] is the intermediate spectral coordinate, in SI units.

  • phi: ndarray

    • Array as double array[ncoord].

  • theta: ndarray

    • Longitude and latitude in the native coordinate system of the projection, in degrees, as double array[ncoord].

  • world: ndarray

    • Array of world coordinates as double array[ncoord][nelem]. For celestial axes, world[][self.lng] and world[][self.lat] are the celestial longitude and latitude, in degrees. For spectral axes, world[][self.spec] is the intermediate spectral coordinate, in SI units.

  • stat: ndarray

    • Status return value for each coordinate as int array[ncoord]. 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.

See also

astropy.wcs.Wcsprm.lat, astropy.wcs.Wcsprm.lng

Definition of the latitude and longitude axes

print_contents()#

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

To get a string of the contents, use repr.

s2p(world, origin)#

Transforms world coordinates to pixel coordinates.

Parameters:
worldndarray

Array of world coordinates, in decimal degrees, as double array[ncoord][nelem].

originint

Specifies the origin of pixel values. The Fortran and FITS standards use an origin of 1. Numpy and C use array indexing with origin at 0.

Returns:
resultdict

Returns a dictionary with the following keys:

  • phi: double array[ncoord]

  • theta: double array[ncoord]

    • Longitude and latitude in the native coordinate system of the projection, in degrees.

  • imgcrd: double array[ncoord][nelem]

    • Array of intermediate world coordinates. For celestial axes, imgcrd[][self.lng] and imgcrd[][self.lat] are the projected x-, and y-coordinates, in pseudo “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: double array[ncoord][nelem]

    • Array of pixel coordinates. Pixel coordinates are zero-based.

  • stat: int array[ncoord]

    • 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.

See also

astropy.wcs.Wcsprm.lat, astropy.wcs.Wcsprm.lng

Definition of the latitude and longitude axes

set()#

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

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

Some attributes that are based on other attributes (such as lattyp on ctype) may not be correct until after set is called.

set strips off trailing blanks in all string members.

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.

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_ps(ps)#

Sets PSi_ma keywords for each i and m.

Parameters:
pssequence of tuple

The input must be a sequence of tuples of the form (i, m, value):

  • i: int. Axis number, as in PSi_ma, (i.e. 1-relative)

  • m: int. Parameter number, as in PSi_ma, (i.e. 0-relative)

  • value: string. Parameter value.

set_pv(pv)#

Sets PVi_ma keywords for each i and m.

Parameters:
pvlist of tuple

The input must be a sequence of tuples of the form (i, m, value):

  • i: int. Axis number, as in PVi_ma, (i.e. 1-relative)

  • m: int. Parameter number, as in PVi_ma, (i.e. 0-relative)

  • value: float. Parameter value.

spcfix() int#

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

Returns:
successint

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

sptr(ctype, i=-1)#

Translates the spectral axis in a WCS object.

For example, a FREQ axis may be translated into ZOPT-F2W and vice versa.

Parameters:
ctypestr

Required spectral CTYPEia, maximum of 8 characters. The first four characters are required to be given and are never modified. The remaining four, the algorithm code, are completely determined by, and must be consistent with, the first four characters. 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.

iint

Index of the spectral axis (0-relative). 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).

sub(axes)#

Extracts the coordinate description for a subimage from a WCS object.

The world coordinate system of the subimage must be separable in the sense that the world coordinates at any point in the subimage must depend only on the pixel coordinates of the axes extracted. In practice, this means that the PCi_ja matrix of the original image must not contain non-zero off-diagonal terms that associate any of the subimage axes with any of the non-subimage axes.

sub can also add axes to a wcsprm object. The new axes will be created using the defaults set by the Wcsprm constructor which produce a simple, unnamed, linear axis with world coordinates equal to the pixel coordinate. These default values can be changed before invoking set.

Parameters:
axesint or a sequence.
  • If an int, include the first N axes in their original order.

  • If a sequence, may contain a combination of image axis numbers (1-relative) or special axis identifiers (see below). Order is significant; axes[0] is the axis number of the input image that corresponds to the first axis in the subimage, etc. Use an axis number of 0 to create a new axis using the defaults.

  • If 0, [] or None, do a deep copy.

Coordinate axes types may be specified using either strings or special integer constants. The available types are:

  • 'longitude' / WCSSUB_LONGITUDE: Celestial longitude

  • 'latitude' / WCSSUB_LATITUDE: Celestial latitude

  • 'cubeface' / WCSSUB_CUBEFACE: Quadcube CUBEFACE axis

  • 'spectral' / WCSSUB_SPECTRAL: Spectral axis

  • 'stokes' / WCSSUB_STOKES: Stokes axis

  • 'temporal' / WCSSUB_TIME: Time axis (requires WCSLIB version 7.8 or greater)

  • 'celestial' / WCSSUB_CELESTIAL: An alias for the combination of 'longitude', 'latitude' and 'cubeface'.

Returns:
new_wcsWCS object
Raises:
MemoryError

Memory allocation failed.

InvalidSubimageSpecificationError

Invalid subimage specification (no spectral axis).

NonseparableSubimageCoordinateSystemError

Non-separable subimage coordinate system.

Notes

Combinations of subimage axes of particular types may be extracted in the same order as they occur in the input image by combining the integer constants with the ‘binary or’ (|) operator. For example:

wcs.sub([WCSSUB_LONGITUDE | WCSSUB_LATITUDE | WCSSUB_SPECTRAL])

would extract the longitude, latitude, and spectral axes in the same order as the input image. If one of each were present, the resulting object would have three dimensions.

For convenience, WCSSUB_CELESTIAL is defined as the combination WCSSUB_LONGITUDE | WCSSUB_LATITUDE | WCSSUB_CUBEFACE.

The codes may also be negated to extract all but the types specified, for example:

wcs.sub([
  WCSSUB_LONGITUDE,
  WCSSUB_LATITUDE,
  WCSSUB_CUBEFACE,
  -(WCSSUB_SPECTRAL | WCSSUB_STOKES)])

The last of these specifies all axis types other than spectral or Stokes. Extraction is done in the order specified by axes, i.e. a longitude axis (if present) would be extracted first (via axes[0]) and not subsequently (via axes[3]). Likewise for the latitude and cubeface axes in this example.

The number of dimensions in the returned object may be less than or greater than the length of axes. However, it will never exceed the number of axes in the input image.

to_header(relax=False)#

to_header translates a WCS object into a FITS header.

The details of the header depends on context:

  • If the colnum member is non-zero then a binary table image array header will be produced.

  • Otherwise, if the colax member is set non-zero then a pixel list header will be produced.

  • Otherwise, a primary image or image extension header will be produced.

The output header will almost certainly differ from the input in a number of respects:

  1. The output header only contains WCS-related keywords. In particular, it does not contain syntactically-required keywords such as SIMPLE, NAXIS, BITPIX, or END.

  2. Deprecated (e.g. CROTAn) or non-standard usage will be translated to standard (this is partially dependent on whether fix was applied).

  3. Quantities will be converted to the units used internally, basically SI with the addition of degrees.

  4. Floating-point quantities may be given to a different decimal precision.

  5. Elements of the PCi_j matrix will be written if and only if they differ from the unit matrix. Thus, if the matrix is unity then no elements will be written.

  6. Additional keywords such as WCSAXES, CUNITia, LONPOLEa and LATPOLEa may appear.

  7. The original keycomments will be lost, although to_header tries hard to write meaningful comments.

  8. Keyword order may be changed.

Keywords can be translated between the image array, binary table, and pixel lists forms by manipulating the colnum or colax members of the WCS object.

Parameters:
relaxbool or int

Degree of permissiveness:

  • False: Recognize only FITS keywords defined by the published WCS standard.

  • True: Admit all recognized informal extensions of the WCS standard.

  • int: a bit field selecting specific extensions to write. See Header-writing relaxation constants for details.

Returns:
headerstr

Raw FITS header as a string.

unitfix(translate_units='')#

Translates non-standard CUNITia keyvalues.

For example, DEG -> deg, also stripping off unnecessary whitespace.

Parameters:
translate_unitsstr, optional

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:
successint

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