ChebyshevField#
- final class lsst.images.fields.ChebyshevField(bounds: Bounds, coefficients: ndarray, *, unit: UnitBase | None = None)#
Bases:
BaseFieldA 2-d Chebyshev polynomial over a rectangular region.
Parameters#
- bounds
The region where this field can be evaluated. The
bboxof this region is grown by half a pixel on all sides and then used to remap coordinates to[-1, 1]x[-1, 1], which is the natural domain of a 2-d Chebyshev polynomial.- coefficients
Coefficients for the 2-d Chebyshev polynomial of the first kind, as a 2-d matrix in which element
[p, q]corresponds to the coefficient ofT_p(y) T_q(x). Will be set to read-only in place.- unit
Units of the field.
Attributes Summary
The region over which this field can be evaluated (
Bounds).Coefficients for the 2-d Chebyshev polynomial of the first kind, as a 2-d matrix in which element
[p, q]corresponds to the coefficient ofT_p(y) T_q(x).Maximum polynomial order in either dimension (
int).The units of the field (
astropy.units.UnitBaseorNone).Maximum polynomial order in the column dimension (
int).Maximum polynomial order in the row dimension (
int).Methods Summary
__call__(*, x, y[, quantity])Call self as a function.
deserialize(model, archive)Deserialize the Chebyshev field from an input archive.
evaluate(*, x, y, quantity)Evaluate at non-gridded points.
fit(bounds, data[, order, weight, y_order, ...])Fit a Chebyshev field to data points using linear least squares.
from_legacy(legacy[, unit])Convert from a legacy
lsst.afw.math.ChebyshevBoundedField.from_legacy_background(legacy_background[, unit])Convert from a legacy
lsst.afw.math.BackgroundMIinstance.multiply_constant(factor)Multiply by a constant, returning a new field of the same type.
render([bbox, dtype])Create an image realization of the field.
serialize(archive)Serialize the Chebyshev field to an output archive.
Convert to a legacy
lsst.afw.math.ChebyshevBoundedField.Attributes Documentation
- bounds#
- coefficients#
Coefficients for the 2-d Chebyshev polynomial of the first kind, as a 2-d matrix in which element
[p, q]corresponds to the coefficient ofT_p(y) T_q(x).
- order#
Maximum polynomial order in either dimension (
int).
- unit#
- x_order#
Maximum polynomial order in the column dimension (
int).
- y_order#
Maximum polynomial order in the row dimension (
int).
Methods Documentation
- __call__(*, x: ndarray, y: ndarray, quantity: bool = False) ndarray | Quantity#
Call self as a function.
- static deserialize(model: ChebyshevFieldSerializationModel, archive: InputArchive[Any]) ChebyshevField#
Deserialize the Chebyshev field from an input archive.
- evaluate(*, x: ndarray, y: ndarray, quantity: bool) ndarray | Quantity#
Evaluate at non-gridded points.
Parameters#
- x
X coordinates to evaluate at.
- y
Y coordinates to evaluate at; must be broadcast-compatible with
x.- quantity
If
True, return anastropy.units.Quantityinstead of anumpy.ndarray. IfunitisNone, the returned object will be a dimensionlessQuantity.
- static fit(bounds: Bounds, data: ndarray | Quantity, order: int | None = None, *, y: ndarray, x: ndarray, weight: ndarray | None = None, y_order: int | None = None, x_order: int | None = None, triangular: bool = True, unit: UnitBase | None = None) ChebyshevField#
Fit a Chebyshev field to data points using linear least squares.
Parameters#
- data
Data points to fit. If this is an
astropy.units.Quantity, this sets the units of the field and theunitargument cannot also be provided.- order
Maximum order for the Chebyshev polynomial in both dimensions.
- y
Y coordinates of the data points. Must have either the same shape as
data(providing the coordinates for all points directly), or be a 1-d array with the same size asdata.shape[0](whendatais a 2-d image andyprovides the coordinates of the rows).- x
X coordinates of the data points. Must have either the same shape as
data(providing the coordinates for all points directly), or be a 1-d array with the same size asdata.shape[1](whendatais a 2-d image andxprovides the coordinates of the columns).- weight
Weights to apply to the data points. Must have the same shape as
data.- y_order
Maximum order for the Chebyshev polynomial in
y. Requiresx_orderto also be provided. Incompatible withorder.- x_order
Maximum order for the Chebyshev polynomial in
x. Requiresy_orderto also be provided. Incompatible withorder.- triangular
If
True, only fit for coefficients ofT_p(y) T_q(x)wherep + q <= max(y_order, x_order).- unit
Units of the returned field.
- static from_legacy(legacy: LegacyChebyshevBoundedField, unit: astropy.units.UnitBase | None = None) ChebyshevField#
Convert from a legacy
lsst.afw.math.ChebyshevBoundedField.
- static from_legacy_background(legacy_background: LegacyBackground, unit: astropy.units.UnitBase | None = None) ChebyshevField#
Convert from a legacy
lsst.afw.math.BackgroundMIinstance.
- multiply_constant(factor: float | Quantity | UnitBase) ChebyshevField#
Multiply by a constant, returning a new field of the same type.
Parameters#
- factor
Factor to multiply by. When this has units, those should multiply
self.unitor set the units of the returned field ifself.unit is None.
- render(bbox: Box | None = None, *, dtype: type[Any] | dtype[Any] | _SupportsDType[dtype[Any]] | tuple[Any, Any] | list[Any] | _DTypeDict | str | None = None) Image#
Create an image realization of the field.
Parameters#
- bbox
Bounding box of the image. If not provided,
self.bounds.bboxwill be used.- dtype
Pixel data type for the returned image.
- serialize(archive: OutputArchive[Any]) ChebyshevFieldSerializationModel#
Serialize the Chebyshev field to an output archive.
- to_legacy() LegacyChebyshevBoundedField#
Convert to a legacy
lsst.afw.math.ChebyshevBoundedField.