Image#

final class lsst.images.Image(array_or_fill: ndarray | int | float = 0, /, *, bbox: Box | None = None, start: Sequence[int] | None = None, shape: Sequence[int] | None = None, dtype: type[Any] | dtype[Any] | _SupportsDType[dtype[Any]] | tuple[Any, Any] | list[Any] | _DTypeDict | str | None = None, unit: UnitBase | None = None, projection: Projection[Any] | None = None, obs_info: ObservationInfo | None = None, metadata: dict[str, MetadataValue] | None = None)#

Bases: GeneralizedImage

A 2-d array that may be augmented with units and a nonzero origin.

Parameters#

array_or_fill

Array or fill value for the image. If a fill value, bbox or shape must be provided.

bbox

Bounding box for the image.

start

Logical coordinates of the first pixel in the array, ordered y, x (unless an XY instance is passed). Ignored if bbox is provided. Defaults to zeros.

shape

Leading dimensions of the array, ordered y, x (unless an XY instance is passed). Only needed if array_or_fill is not an array and bbox is not provided. Like the bbox, this does not include the last dimension of the array.

dtype

Pixel data type override.

unit

Units for the image’s pixel values.

projection

Projection that maps the pixel grid to the sky.

obs_info

General information about the associated observation in standardized form.

metadata

Arbitrary flexible metadata to associate with the image.

Notes#

Indexing the array attribute of an Image does not take into account its start offset, but accessing a subimage by indexing an Image with a Box does, and the bbox of the subimage is set to match its location within the original image.

Indexed assignment to a subimage requires consistency between the coordinate systems and units of both operands, but it will automatically select a subimage of the right-hand side and convert compatible units when possible. In other words:

a[box] = b

is a shortcut for

a[box].quantity = b[box].quantity

An ellipsis (...) can be used instead of a Box to assign to the full image.

Attributes Summary

absolute

A proxy object for slicing a generalized image using absolute pixel coordinates.

array

The low-level array (numpy.ndarray).

astropy_wcs

An Astropy WCS for this image's pixel array.

bbox

Bounding box for the image (Box).

butler_dataset

The butler dataset reference for this image (lsst.daf.butler.SerializedDatasetRef | None).

butler_provenance

The butler inputs and ID of the task quantum that produced this dataset (lsst.daf.butler.DatasetProvenance | None)

fits_wcs

An Astropy FITS WCS for this image's pixel array.

local

A proxy object for slicing a generalized image using "local" or "array" pixel coordinates.

metadata

Arbitrary flexible metadata associated with the image (dict).

obs_info

General information about the associated observation in standard form.

projection

The projection that maps this image's pixel grid to the sky (Projection | None).

quantity

The low-level array with units (astropy.units.Quantity).

unit

Units for the image's pixel values (astropy.units.Unit or None).

Methods Summary

copy()

Deep-copy the image and metadata.

deserialize(model, archive, *[, bbox, ...])

Deserialize an image from an input archive.

from_legacy(legacy[, unit])

Convert from an lsst.afw.image.Image instance.

read_fits(url, *[, bbox])

Read an image from a FITS file.

read_legacy(uri, *[, preserve_quantization, ...])

Read a FITS file written by lsst.afw.image.Image.writeFits.

serialize(archive, *[, update_header, ...])

Serialize the image to an output archive.

to_legacy(*[, copy])

Convert to an lsst.afw.image.Image instance.

view(*[, unit, projection, start, obs_info])

Make a view of the image, with optional updates.

write_fits(filename, *[, compression, ...])

Write the image to a FITS file.

Attributes Documentation

absolute#

A proxy object for slicing a generalized image using absolute pixel coordinates.

Notes#

In this convention, the first row and column of the pixel grid is bbox.start. A subimage and its parent image share the same absolute pixel coordinate system, and most lsst.images types (e.g. Box, Projection, PointSpreadFunction) operate exclusively in this system.

Note that astropy.wcs and numpy.ndarray are not aware of the bbox.start offset that defines tihs coordinates system; use local slicing for indices obtained from those.

See Also#

lsst.images.BoxSliceFactory lsst.images.IntervalSliceFactory

array#

The low-level array (numpy.ndarray).

Assigning to this attribute modifies the existing array in place; the bounding box and underlying data pointer are never changed.

astropy_wcs#

An Astropy WCS for this image’s pixel array.

Notes#

As expected for Astropy WCS objects, this defines pixel coordinates such that the first row and column in any associated arrays are (0, 0), not bbox.start, as is the case for projection.

This object satisfies the astropy.wcs.wcsapi.BaseHighLevelWCS and astropy.wcs.wcsapi.BaseLowLevelWCS interfaces, but it is not an astropy.wcs.WCS (use fits_wcs for that).

bbox#

Bounding box for the image (Box).

butler_dataset#

The butler dataset reference for this image (lsst.daf.butler.SerializedDatasetRef | None).

butler_provenance#

The butler inputs and ID of the task quantum that produced this dataset (lsst.daf.butler.DatasetProvenance | None)

fits_wcs#

An Astropy FITS WCS for this image’s pixel array.

Notes#

As expected for Astropy WCS objects, this defines pixel coordinates such that the first row and column in any associated arrays are (0, 0), not bbox.start, as is the case for projection.

This may be an approximation or absent if projection is not naturally representable as a FITS WCS.

local#

A proxy object for slicing a generalized image using “local” or “array” pixel coordinates.

Notes#

In this convention, the first row and column of the pixel grid is always at (0, 0). This is also the convention used by astropy.wcs objects. When a subimage is created from a parent image, its “local” coordinate system is offset from the coordinate systems of the parent image.

Note that most lsst.images types (e.g. Box, Projection, PointSpreadFunction) operate instead in “absolute” coordinates, which is shared by subimage and their parents.

See Also#

lsst.images.BoxSliceFactory lsst.images.IntervalSliceFactory

metadata#

Arbitrary flexible metadata associated with the image (dict).

Notes#

Metadata is shared with subimages and other views. It can be disconnected by reassigning to a copy explicitly:

image.metadata = image.metadata.copy()

obs_info#

General information about the associated observation in standard form. (ObservationInfo | None).

projection#

The projection that maps this image’s pixel grid to the sky (Projection | None).

Notes#

The pixel coordinates used by this projection account for the bounding box start; they are not just array indices.

quantity#

The low-level array with units (astropy.units.Quantity).

Assigning to this attribute modifies the existing array in place; the bounding box and underlying data pointer are never changed.

unit#

Units for the image’s pixel values (astropy.units.Unit or None).

Methods Documentation

copy() Image#

Deep-copy the image and metadata.

Attached immutable objects (like Projection instances) are not copied.

static deserialize(model: ~lsst.images._image.ImageSerializationModel[Any], archive: ~lsst.images.serialization._input_archive.InputArchive[~typing.Any], *, bbox: ~lsst.images._geom.Box | None = None, strip_header: ~collections.abc.Callable[[~astropy.io.fits.header.Header], None] = <function no_header_updates>) Image#

Deserialize an image from an input archive.

Parameters#

model

A Pydantic model representation of the image, holding references to data stored in the archive.

archive

Archive to read from.

bbox

Bounding box of a subimage to read instead.

strip_header

A callable that strips out any FITS header cards added by the update_header argument in the corresponding call to serialize.

static from_legacy(legacy: Any, unit: UnitBase | None = None) Image#

Convert from an lsst.afw.image.Image instance.

Parameters#

legacy

An lsst.afw.image.Image instance that will share pixel data with

the returned object.

unit

Units of the image.

static read_fits(url: str | ParseResult | ResourcePath | Path, *, bbox: Box | None = None) Image#

Read an image from a FITS file.

Parameters#

url

URL of the file to read; may be any type supported by lsst.resources.ResourcePath.

bbox

Bounding box of a subimage to read instead.

static read_legacy(uri: str | ParseResult | ResourcePath | Path, *, preserve_quantization: bool = False, ext: str | int = 1, fits_wcs_frame: Frame | None = None) Image#

Read a FITS file written by lsst.afw.image.Image.writeFits.

Parameters#

uri

URI or file name.

preserve_quantization

If True, ensure that writing the image back out again will exactly preserve quantization-compressed pixel values. This causes the arrays to be marked as read-only and stores the original binary table data for those planes in memory. If the Image is copied, the precompressed pixel values are not transferred to the copy.

ext

Name or index of the FITS HDU to read.

fits_wcs_frame

If not None and the HDU containing the image has a FITS WCS, attach a Projection to the returned image by converting that WCS.

serialize(archive: ~lsst.images.serialization._output_archive.OutputArchive, *, update_header: ~collections.abc.Callable[[~astropy.io.fits.header.Header], None] = <function no_header_updates>, save_projection: bool = True, save_obs_info: bool = True, add_offset_wcs: str | None = 'A') ImageSerializationModel[TypeVar]#

Serialize the image to an output archive.

Parameters#

archive

Archive to write to.

update_header

A callback that will be given the FITS header for the HDU containing this image in order to add keys to it. This callback may be provided but will not be called if the output format is not FITS.

save_projection

If True, save the Projection attached to the image, if there is one. This does not affect whether a FITS WCS corresponding to the projection is written (it always is, if available, and if add_offset_wcs is not " ").

save_obs_info

If True, save the ObservationInfo attached to the image, if there is one.

add_offset_wcs

A FITS WCS single-character suffix to use when adding a linear WCS that maps the FITS array to the logical pixel coordinates defined by bbox.start. Set to None to not write this WCS. If this is set to " ", it will prevent the Projection from being saved as a FITS WCS.

to_legacy(*, copy: bool | None = None) Any#

Convert to an lsst.afw.image.Image instance.

Parameters#

copy

If True, always copy the pixel data. If False, return a view, and raise TypeError if the pixel data is read-only (this is not supported by afw). If None, onyl if the pixel data is read-only.

view(*, unit: UnitBase | None | ellipsis = Ellipsis, projection: Projection | None | ellipsis = Ellipsis, start: Sequence[int] | ellipsis = Ellipsis, obs_info: ObservationInfo | None | ellipsis = Ellipsis) Image#

Make a view of the image, with optional updates.

write_fits(filename: str, *, compression: ~lsst.images.fits._common.FitsCompressionOptions | None = FitsCompressionOptions(algorithm=<FitsCompressionAlgorithm.GZIP_2: 'GZIP_2'>, tile_shape=None, quantization=None), compression_seed: int | None = None) None#

Write the image to a FITS file.

Parameters#

filename

Name of the file to write to. Must be a local file.

compression

Compression options.

compression_seed

A FITS tile compression seed to use whenever the configured compression seed is None or (for backwards compatibility) 0.