Bases: object
Wcsprm is a direct wrapper around wcslib. It provides access to the core WCS transformations that it supports.
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.
Parameters: | header : An astropy.io.fits.Header, string, or None.
key : str, optional
relax : bool or int, optional
naxis : int, optional
keysel : sequence of flag bits, optional
colsel : sequence of int
|
---|---|
Raises: | MemoryError
ValueError
KeyError
|
Attributes Summary
alt | str Character code for alternate coordinate descriptions. |
axis_types | int array[naxis] An array of four-digit type codes for each axis. |
cd | double array[naxis][naxis] The CDi_ja linear transformation |
cdelt | double array[naxis] Coordinate increments (CDELTia) for each |
cel_offset | boolean Is there an offset? |
cname | list of strings A list of the coordinate axis names, from |
colax | int array[naxis] An array recording the column numbers for each |
colnum | int Column of FITS binary table associated with this WCS. |
crder | double array[naxis] The random error in each coordinate axis, |
crota | double array[naxis] CROTAia keyvalues for each coordinate |
crpix | double array[naxis] Coordinate reference pixels (CRPIXja) for |
crval | double array[naxis] Coordinate reference values (CRVALia) for |
csyer | double array[naxis] The systematic error in the coordinate value |
ctype | list of strings[naxis] List of CTYPEia keyvalues. |
cubeface | int Index into the pixcrd (pixel coordinate) array for the |
cunit | list of astropy.UnitBase[naxis] List of CUNITia keyvalues as |
dateavg | string Representative mid-point of the date of observation. |
dateobs | string Start of the date of observation. |
equinox | double The equinox associated with dynamical equatorial or |
imgpix_matrix | double array[2][2] (read-only) Inverse of the CDELT or PC |
lat | int (read-only) The index into the world coord array containing |
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 |
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. |
mjdobs | double Modified Julian Date corresponding to DATE-OBS. |
name | string The name given to the coordinate representation |
naxis | int (read-only) The number of axes (pixel and coordinate). |
obsgeo | double array[3] Location of the observer in a standard terrestrial |
pc | double array[naxis][naxis] The PCi_ja (pixel coordinate) |
phi0 | double The native latitude of the fiducial point. |
piximg_matrix | double array[2][2] (read-only) Matrix containing the product of |
radesys | string The equatorial or ecliptic coordinate system type, |
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. |
theta0 | double The native longitude of the fiducial point. |
velangl | double Velocity angle. |
velosys | double Relative radial velocity. |
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]) | 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(() -> double array[naxis]) | Coordinate increments (CDELTia) for each coord axis. |
get_pc(() -> double array[naxis][naxis]) | Returns the PC matrix in read-only form. |
get_ps(() -> list of tuples) | Returns PSi_ma keywords for each i and m. |
get_pv(() -> list of tuples) | Returns PVi_ma keywords for each i and m. |
has_cd(() -> bool) | Returns True if CDi_ja is present. |
has_cdi_ja(() -> bool) | Alias for has_cd. |
has_crota(() -> bool) | Returns True if CROTAia is present. |
has_crotaia(() -> bool) | Alias for has_crota. |
has_pc(() -> bool) | Returns True if PCi_ja is present. |
has_pci_ja(() -> bool) | Alias for has_pc. |
is_unity(() -> bool) | 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(() -> int) | 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
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.
int array[naxis] An array of four-digit type codes for each axis.
CTYPEia in "4-3" form with unrecognized algorithm code will have its type set to -1 and generate an error.
double array[naxis][naxis] The CDi_ja linear transformation matrix.
For historical compatibility, three alternate specifications of the linear transforations 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.
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.
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.
list of strings A list of the coordinate axis names, from CNAMEia.
int array[naxis] An array recording the column numbers for each axis in a pixel list.
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.
double array[naxis] The random error in each coordinate axis, CRDERia.
An undefined value is represented by NaN.
double array[naxis] CROTAia keyvalues for each coordinate axis.
For historical compatibility, three alternate specifications of the linear transforations 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.
double array[naxis] Coordinate reference pixels (CRPIXja) for each pixel axis.
double array[naxis] Coordinate reference values (CRVALia) for each coordinate axis.
double array[naxis] The systematic error in the coordinate value axes, CSYERia.
An undefined value is represented by NaN.
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.
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 5Faces 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.
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.
string Representative mid-point of the date of observation.
In ISO format, yyyy-mm-ddThh:mm:ss.
See also
string Start of the date of observation.
In ISO format, yyyy-mm-ddThh:mm:ss.
See also
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.
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.
int (read-only) The index into the world coord array containing latitude values.
double The native latitude of the celestial pole, LATPOLEa (deg).
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.
int (read-only) The index into the world coord array containing longitude values.
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.
double The native longitude of the celestial pole.
LONPOLEa (deg).
double Modified Julian Date corresponding to DATE-AVG.
(MJD = JD - 2400000.5).
An undefined value is represented by NaN.
See also
double Modified Julian Date corresponding to DATE-OBS.
(MJD = JD - 2400000.5).
An undefined value is represented by NaN.
See also
string The name given to the coordinate representation WCSNAMEa.
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:
- NAXIS
- WCSAXESa
- 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.
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.
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 transforations 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.
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.
See also
double array[2][2] (read-only) Matrix containing the product of the CDELTia diagonal matrix and the PCi_ja matrix.
string The equatorial or ecliptic coordinate system type, RADESYSa.
double Rest frequency (Hz) from RESTFRQa.
An undefined value is represented by NaN.
double Rest wavelength (m) from RESTWAVa.
An undefined value is represented by NaN.
int (read-only) The index containing the spectral axis values.
string Spectral reference frame (standard of rest), SPECSYSa.
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.
string Spectral reference frame for redshift.
The spectral reference frame (standard of rest) in which the redshift was measured, SSYSSRCa.
list of Tabprm Tabular coordinate objects.
A list of tabular coordinate objects associated with this WCS.
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.
See also
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.
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.
double The redshift, ZSOURCEa, of the source.
An undefined value is represented by NaN.
Methods Documentation
Enable/disable bounds checking.
Parameters: | pix2world : bool, optional world2pix : bool, optional |
---|
Notes
Note that by default (without calling bounds_check) strict bounds checking is enabled.
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: | success : int
|
---|
Translates AIPS-convention celestial projection types, -NCP and -GLS.
Returns: | success : int
|
---|
Compare two Wcsprm objects for equality.
Parameters: | other : Wcsprm
cmp : int, optional
|
---|---|
Returns: | equal : bool |
Fixes WCS keyvalues for malformed cylindrical projections.
Returns: | success : int
|
---|
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: | success : int
|
---|
Applies all of the corrections handled separately by datfix, unitfix, celfix, spcfix, cylfix and cdfix.
Parameters: | translate_units : str, optional
naxis : int array[naxis], optional
|
---|---|
Returns: | status : dict |
Coordinate increments (CDELTia) for each coord axis.
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.
Returns the PC matrix in read-only form. 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.
Returns PSi_ma keywords for each i and m.
Returns: | ps : list of tuples
|
---|
See also
Returns PVi_ma keywords for each i and m.
Returns: | Returned as a list of tuples of the form (i, m, value):
|
---|
See also
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.
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
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
Returns True if PCi_ja is present. PCi_ja is the recommended way to specify the linear transformation matrix.
See also
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
mixcel : int
vspan : pair of floats
vstep : float
viter : int
world : double array[naxis]
pixcrd : double array[naxis].
origin : int
|
---|---|
Returns: | result : dict
|
Raises: | MemoryError
SingularMatrixError
InconsistentAxisTypesError
ValueError
InvalidTransformError
InvalidTransformError
InvalidCoordinateError
NoSolutionError
|
See also
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.
Converts pixel to world coordinates.
Parameters: | pixcrd : double array[ncoord][nelem]
origin : int
|
---|---|
Returns: | result : dict
|
Raises: | MemoryError
SingularMatrixError
InconsistentAxisTypesError
ValueError
ValueError
InvalidTransformError
InvalidTransformError
|
See also
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.
Transforms world coordinates to pixel coordinates.
Parameters: | world : double array[ncoord][nelem]
origin : int
|
---|---|
Returns: | result : dict
|
Raises: | MemoryError
SingularMatrixError
InconsistentAxisTypesError
ValueError
InvalidTransformError
InvalidTransformError
|
See also
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
SingularMatrixError
InconsistentAxisTypesError
ValueError
InvalidTransformError
InvalidTransformError
|
---|
Sets PSi_ma keywords for each i and m.
Parameters: | ps : sequence of tuples
|
---|
See also
Sets PVi_ma keywords for each i and m.
Parameters: | pv : list of tuples
|
---|
See also
Translates AIPS-convention spectral coordinate types. {FREQ, VELO, FELO}-{OBS, HEL, LSR} (e.g. FREQ-LSR, VELO-OBS, FELO-HEL)
Returns: | success : int
|
---|
Translates the spectral axis in a WCS object.
For example, a FREQ axis may be translated into ZOPT-F2W and vice versa.
Parameters: | ctype : str
i : int
|
---|---|
Raises: | MemoryError
SingularMatrixError
InconsistentAxisTypesError
ValueError
InvalidTransformError
InvalidTransformError
InvalidSubimageSpecificationError
|
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: | axes : int or a sequence.
|
---|---|
Returns: | new_wcs : WCS object |
Raises: | MemoryError
InvalidSubimageSpecificationError
NonseparableSubimageCoordinateSystem
|
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 translates a WCS object into a FITS header.
The details of the header depends on context:
The output header will almost certainly differ from the input in a number of respects:
- The output header only contains WCS-related keywords. In particular, it does not contain syntactically-required keywords such as SIMPLE, NAXIS, BITPIX, or END.
- Deprecated (e.g. CROTAn) or non-standard usage will be translated to standard (this is partially dependent on whether fix was applied).
- Quantities will be converted to the units used internally, basically SI with the addition of degrees.
- Floating-point quantities may be given to a different decimal precision.
- 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.
- Additional keywords such as WCSAXES, CUNITia, LONPOLEa and LATPOLEa may appear.
- The original keycomments will be lost, although to_header tries hard to write meaningful comments.
- 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: | relax : bool or int
|
---|---|
Returns: | header : str
|
Translates non-standard CUNITia keyvalues.
For example, DEG -> deg, also stripping off unnecessary whitespace.
Parameters: | translate_units : str, optional
|
---|---|
Returns: | success : int
|