GriddedPSFModel#

class photutils.psf.GriddedPSFModel(nddata, *, flux=1.0, x_0=0.0, y_0=0.0, fill_value=0.0)[source]#

Bases: ModelGridPlotMixin, Fittable2DModel

A model for a grid of 2D ePSF models.

The ePSF models are defined at fiducial detector locations and are bilinearly interpolated to calculate an ePSF model at an arbitrary (x, y) detector position. The fiducial detector locations are must form a rectangular grid.

The model has three model parameters: an image intensity scaling factor (flux) which is applied to the input image, and two positional parameters (x_0 and y_0) indicating the location of a feature in the coordinate grid on which the model is evaluated.

When evaluating this model, it cannot be called with x and y arrays that have greater than 2 dimensions.

Parameters:
nddataNDData

A NDData object containing the grid of reference ePSF arrays. The data attribute must contain a 3D ndarray containing a stack of the 2D ePSFs with a shape of (N_psf, ePSF_ny, ePSF_nx). The length of the x and y axes must both be at least 4. All elements of the input image data must be finite. The PSF peak is assumed to be located at the center of the input image. The array must be normalized so that the total flux of a source is 1.0. This means that the sum of the values in the input image PSF over an infinite grid is 1.0. In practice, the sum of the data values in the input image may be less than 1.0 if the input image only covers a finite region of the PSF. These correction factors can be estimated from the ensquared or encircled energy of the PSF based on the size of the input image.

The meta attribute must be dictionary containing the following:

  • 'grid_xypos': A list of the (x, y) grid positions of each reference ePSF. The order of positions should match the first axis of the 3D ndarray of ePSFs. In other words, grid_xypos[i] should be the (x, y) position of the reference ePSF defined in nddata.data[i]. The grid positions must form a rectangular grid.

  • 'oversampling': The integer oversampling factor(s) of the input ePSF images. If oversampling is a scalar then it will be used for both axes. If oversampling has two elements, they must be in (y, x) order.

The meta attribute may contain other properties such as the telescope, instrument, detector, and filter of the ePSF.

fluxfloat, optional

The flux scaling factor for the model. This is the total flux of the source, assuming the input ePSF images are properly normalized.

x_0, y_0float, optional

The (x, y) position of the PSF peak in the image in the output coordinate grid on which the model is evaluated.

fill_valuefloat, optional

The value to use for points outside of the input pixel grid. The default is 0.0.

Methods

read(*args, **kwargs)

Class method to create a GriddedPSFModel instance from a STDPSF FITS file. This method uses stdpsf_reader() with the provided parameters.

Notes

Internally, the grid of ePSFs will be arranged and stored such that it is sorted first by the y reference pixel coordinate and then by the x reference pixel coordinate.

Attributes Summary

bbox_with_units

bounding_box

A tuple of length n_inputs defining the bounding box limits, or raise NotImplementedError for no bounding_box.

bounds

A dict mapping parameter names to their upper and lower bounds as (min, max) tuples or [min, max] lists.

col_fit_deriv

cov_matrix

Fitter should set covariance matrix, if available.

eqcons

List of parameter equality constraints.

fit_deriv

Function (similar to the model's evaluate) to compute the derivatives of the model with respect to its parameters, for use by fitting algorithms.

fittable

fixed

A dict mapping parameter names to their fixed constraint.

flux

has_bounds

Whether the model has any bounds constraints.

has_fixed

Whether the model has any fixed constraints.

has_tied

Whether the model has any tied constraints.

has_user_bounding_box

A flag indicating whether or not a custom bounding_box has been assigned to this model by a user, via assignment to model.bounding_box.

has_user_inverse

A flag indicating whether or not a custom inverse model has been assigned to this model by a user, via assignment to model.inverse.

ineqcons

List of parameter inequality constraints.

input_units

This property is used to indicate what units or sets of units the evaluate method expects, and returns a dictionary mapping inputs to units (or None if any units are accepted).

input_units_allow_dimensionless

Allow dimensionless input (and corresponding output).

input_units_equivalencies

input_units_strict

Enforce strict units on inputs to evaluate.

inputs

inverse

Returns a new Model instance which performs the inverse transform, if an analytic inverse is defined for this model.

linear

meta

A dict-like object to store optional information.

model_constraints

Primarily for informational purposes, these are the types of constraints that constrain model evaluation.

model_set_axis

The index of the model set axis--that is the axis of a parameter array that pertains to which model a parameter value pertains to--as specified when the model was initialized.

n_inputs

n_outputs

n_submodels

Return the number of components in a single model, which is obviously 1.

name

User-provided name for this model instance.

origin

A 1D ndarray (x, y) pixel coordinates within the model's 2D image of the origin of the coordinate system.

outputs

oversampling

The integer oversampling factor(s) of the input ePSF images.

param_names

Names of the parameters that describe models of this type.

param_sets

Return parameters as a pset.

parameter_constraints

Primarily for informational purposes, these are the types of constraints that can be set on a model's parameters.

parameters

A flattened array of all parameter values in all parameter sets.

read

Read and parse a FITS file into a GriddedPSFModel instance.

return_units

This property is used to indicate what units or sets of units the output of evaluate should be in, and returns a dictionary mapping outputs to units (or None if any units are accepted).

separable

A flag indicating whether a model is separable.

standard_broadcasting

stds

Standard deviation of parameters, if covariance matrix is available.

sync_constraints

This is a boolean property that indicates whether or not accessing constraints automatically check the constituent models current values.

tied

A dict mapping parameter names to their tied constraint.

uses_quantity

True if this model has been created with Quantity objects or if there are no parameters.

x_0

y_0

Methods Summary

__call__(*inputs[, model_set_axis, ...])

Evaluate this model using the given input(s) and the parameter values that were specified when the model was instantiated.

coerce_units([input_units, return_units, ...])

Attach units to this (unitless) model.

copy()

Return a copy of this model where only the model parameters are copied.

deepcopy()

Return a deep copy of this model.

evaluate(x, y, flux, x_0, y_0)

Calculate the ePSF model at the input coordinates for the given model parameters.

get_bounding_box([with_bbox])

Return the bounding_box of a model if it exists or None otherwise.

has_inverse()

Returns True if the model has an analytic or user inverse defined.

input_shape(inputs)

Get input shape for bounding_box evaluation.

output_units(**kwargs)

Return a dictionary of output units for this model given a dictionary of fitting inputs and outputs.

plot_grid(*[, ax, vmax_scale, peak_norm, ...])

Plot the grid of ePSF models.

prepare_inputs(*inputs[, model_set_axis, ...])

This method is used in __call__ to ensure that all the inputs to the model can be broadcast into compatible shapes (if one or both of them are input as arrays), particularly if there are more than one parameter sets.

prepare_outputs(broadcasted_shapes, ...)

rename([name, inputs, outputs])

Return a copy of this model with a new name.

render([out, coords])

Evaluate a model at fixed positions, respecting the bounding_box.

set_slice_args(*args)

strip_units_from_tree()

sum_of_implicit_terms(*args, **kwargs)

Evaluate the sum of any implicit model terms on some input variables.

with_units_from_data(**kwargs)

Return an instance of the model which has units for which the parameter values are compatible with the data units specified.

without_units_for_data(**kwargs)

Return an instance of the model for which the parameter values have been converted to the right units for the data, then the units have been stripped away.

Attributes Documentation

bbox_with_units#
bounding_box#

A tuple of length n_inputs defining the bounding box limits, or raise NotImplementedError for no bounding_box.

The default limits are given by a bounding_box property or method defined in the class body of a specific model. If not defined then this property just raises NotImplementedError by default (but may be assigned a custom value by a user). bounding_box can be set manually to an array-like object of shape (model.n_inputs, 2). For further usage, see Efficient Model Rendering with Bounding Boxes

The limits are ordered according to the numpy 'C' indexing convention, and are the reverse of the model input order, e.g. for inputs ('x', 'y', 'z'), bounding_box is defined:

  • for 1D: (x_low, x_high)

  • for 2D: ((y_low, y_high), (x_low, x_high))

  • for 3D: ((z_low, z_high), (y_low, y_high), (x_low, x_high))

Examples

Setting the bounding_box limits for a 1D and 2D model:

>>> from astropy.modeling.models import Gaussian1D, Gaussian2D
>>> model_1d = Gaussian1D()
>>> model_2d = Gaussian2D(x_stddev=1, y_stddev=1)
>>> model_1d.bounding_box = (-5, 5)
>>> model_2d.bounding_box = ((-6, 6), (-5, 5))

Setting the bounding_box limits for a user-defined 3D custom_model:

>>> from astropy.modeling.models import custom_model
>>> def const3d(x, y, z, amp=1):
...    return amp
...
>>> Const3D = custom_model(const3d)
>>> model_3d = Const3D()
>>> model_3d.bounding_box = ((-6, 6), (-5, 5), (-4, 4))

To reset bounding_box to its default limits just delete the user-defined value–this will reset it back to the default defined on the class:

>>> del model_1d.bounding_box

To disable the bounding box entirely (including the default), set bounding_box to None:

>>> model_1d.bounding_box = None
>>> model_1d.bounding_box  
Traceback (most recent call last):
NotImplementedError: No bounding box is defined for this model
(note: the bounding box was explicitly disabled for this model;
use `del model.bounding_box` to restore the default bounding box,
if one is defined for this model).
bounds#

A dict mapping parameter names to their upper and lower bounds as (min, max) tuples or [min, max] lists.

col_fit_deriv = True#
cov_matrix#

Fitter should set covariance matrix, if available.

eqcons#

List of parameter equality constraints.

fit_deriv = None#

Function (similar to the model’s evaluate) to compute the derivatives of the model with respect to its parameters, for use by fitting algorithms. In other words, this computes the Jacobian matrix with respect to the model’s parameters.

fittable = True#
fixed#

A dict mapping parameter names to their fixed constraint.

flux = Parameter('flux', value=1.0)#
has_bounds#

Whether the model has any bounds constraints.

has_fixed#

Whether the model has any fixed constraints.

has_tied#

Whether the model has any tied constraints.

has_user_bounding_box#

A flag indicating whether or not a custom bounding_box has been assigned to this model by a user, via assignment to model.bounding_box.

has_user_inverse#

A flag indicating whether or not a custom inverse model has been assigned to this model by a user, via assignment to model.inverse.

ineqcons#

List of parameter inequality constraints.

input_units#

This property is used to indicate what units or sets of units the evaluate method expects, and returns a dictionary mapping inputs to units (or None if any units are accepted).

Model sub-classes can also use function annotations in evaluate to indicate valid input units, in which case this property should not be overridden since it will return the input units based on the annotations.

input_units_allow_dimensionless#

Allow dimensionless input (and corresponding output). If this is True, input values to evaluate will gain the units specified in input_units. If this is a dictionary then it should map input name to a bool to allow dimensionless numbers for that input. Only has an effect if input_units is defined.

input_units_equivalencies = None#
input_units_strict#

Enforce strict units on inputs to evaluate. If this is set to True, input values to evaluate will be in the exact units specified by input_units. If the input quantities are convertible to input_units, they are converted. If this is a dictionary then it should map input name to a bool to set strict input units for that parameter.

inputs#
inverse#

Returns a new Model instance which performs the inverse transform, if an analytic inverse is defined for this model.

Even on models that don’t have an inverse defined, this property can be set with a manually-defined inverse, such a pre-computed or experimentally determined inverse (often given as a PolynomialModel, but not by requirement).

A custom inverse can be deleted with del model.inverse. In this case the model’s inverse is reset to its default, if a default exists (otherwise the default is to raise NotImplementedError).

Note to authors of Model subclasses: To define an inverse for a model simply override this property to return the appropriate model representing the inverse. The machinery that will make the inverse manually-overridable is added automatically by the base class.

linear = False#
meta = None#

A dict-like object to store optional information.

model_constraints = ('eqcons', 'ineqcons')#

Primarily for informational purposes, these are the types of constraints that constrain model evaluation.

model_set_axis#

The index of the model set axis–that is the axis of a parameter array that pertains to which model a parameter value pertains to–as specified when the model was initialized.

See the documentation on Model Sets for more details.

n_inputs = 2#
n_outputs = 1#
n_submodels#

Return the number of components in a single model, which is obviously 1.

name#

User-provided name for this model instance.

origin#

A 1D ndarray (x, y) pixel coordinates within the model’s 2D image of the origin of the coordinate system.

outputs#
oversampling#

The integer oversampling factor(s) of the input ePSF images.

If oversampling is a scalar then it will be used for both axes. If oversampling has two elements, they must be in (y, x) order.

param_names = ('flux', 'x_0', 'y_0')#

Names of the parameters that describe models of this type.

The parameters in this tuple are in the same order they should be passed in when initializing a model of a specific type. Some types of models, such as polynomial models, have a different number of parameters depending on some other property of the model, such as the degree.

When defining a custom model class the value of this attribute is automatically set by the Parameter attributes defined in the class body.

param_sets#

Return parameters as a pset.

This is a list with one item per parameter set, which is an array of that parameter’s values across all parameter sets, with the last axis associated with the parameter set.

parameter_constraints = ('fixed', 'tied', 'bounds')#

Primarily for informational purposes, these are the types of constraints that can be set on a model’s parameters.

parameters#

A flattened array of all parameter values in all parameter sets.

Fittable parameters maintain this list and fitters modify it.

read#

Read and parse a FITS file into a GriddedPSFModel instance.

This class enables the astropy unified I/O layer for GriddedPSFModel. This allows easily reading a file in different supported data formats using syntax such as:

>>> from photutils.psf import GriddedPSFModel
>>> psf_model = GriddedPSFModel.read('filename.fits', format=format)

Get help on the available readers for GriddedPSFModel using the help() method:

>>> # Get help reading Table and list supported formats
>>> GriddedPSFModel.read.help()

>>> # Get detailed help on the STSPSF FITS reader
>>> GriddedPSFModel.read.help('stdpsf')

>>> # Get detailed help on the WebbPSF FITS reader
>>> GriddedPSFModel.read.help('webbpsf')

>>> # Print list of available formats
>>> GriddedPSFModel.read.list_formats()
Parameters:
instanceobject

Descriptor calling instance or None if no instance.

clstype

Descriptor calling class (either owner class or instance class).

return_units#

This property is used to indicate what units or sets of units the output of evaluate should be in, and returns a dictionary mapping outputs to units (or None if any units are accepted).

Model sub-classes can also use function annotations in evaluate to indicate valid output units, in which case this property should not be overridden since it will return the return units based on the annotations.

separable#

A flag indicating whether a model is separable.

standard_broadcasting = True#
stds#

Standard deviation of parameters, if covariance matrix is available.

sync_constraints#

This is a boolean property that indicates whether or not accessing constraints automatically check the constituent models current values. It defaults to True on creation of a model, but for fitting purposes it should be set to False for performance reasons.

tied#

A dict mapping parameter names to their tied constraint.

uses_quantity#

True if this model has been created with Quantity objects or if there are no parameters.

This can be used to determine if this model should be evaluated with Quantity or regular floats.

x_0 = Parameter('x_0', value=0.0)#
y_0 = Parameter('y_0', value=0.0)#

Methods Documentation

__call__(*inputs, model_set_axis=None, with_bounding_box=False, fill_value=nan, equivalencies=None, inputs_map=None, **new_inputs)#

Evaluate this model using the given input(s) and the parameter values that were specified when the model was instantiated.

coerce_units(input_units=None, return_units=None, input_units_equivalencies=None, input_units_allow_dimensionless=False)#

Attach units to this (unitless) model.

Parameters:
input_unitsdict or tuple, optional

Input units to attach. If dict, each key is the name of a model input, and the value is the unit to attach. If tuple, the elements are units to attach in order corresponding to inputs.

return_unitsdict or tuple, optional

Output units to attach. If dict, each key is the name of a model output, and the value is the unit to attach. If tuple, the elements are units to attach in order corresponding to outputs.

input_units_equivalenciesdict, optional

Default equivalencies to apply to input values. If set, this should be a dictionary where each key is a string that corresponds to one of the model inputs.

input_units_allow_dimensionlessbool or dict, optional

Allow dimensionless input. If this is True, input values to evaluate will gain the units specified in input_units. If this is a dictionary then it should map input name to a bool to allow dimensionless numbers for that input.

Returns:
CompoundModel

A CompoundModel composed of the current model plus UnitsMapping model(s) that attach the units.

Raises:
ValueError

If the current model already has units.

Examples

Wrapping a unitless model to require and convert units:

>>> from astropy.modeling.models import Polynomial1D
>>> from astropy import units as u
>>> poly = Polynomial1D(1, c0=1, c1=2)
>>> model = poly.coerce_units((u.m,), (u.s,))
>>> model(u.Quantity(10, u.m))  
<Quantity 21. s>
>>> model(u.Quantity(1000, u.cm))  
<Quantity 21. s>
>>> model(u.Quantity(10, u.cm))  
<Quantity 1.2 s>

Wrapping a unitless model but still permitting unitless input:

>>> from astropy.modeling.models import Polynomial1D
>>> from astropy import units as u
>>> poly = Polynomial1D(1, c0=1, c1=2)
>>> model = poly.coerce_units((u.m,), (u.s,), input_units_allow_dimensionless=True)
>>> model(u.Quantity(10, u.m))  
<Quantity 21. s>
>>> model(10)  
<Quantity 21. s>
copy()[source]#

Return a copy of this model where only the model parameters are copied.

All other copied model attributes are references to the original model. This prevents copying the ePSF grid data, which may contain a large array.

This method is useful if one is interested in only changing the model parameters in a model copy. It is used in the PSF photometry classes during model fitting.

Use the deepcopy method if you want to copy all of the model attributes, including the ePSF grid data.

Returns:
resultGriddedPSFModel

A copy of this model with only the model parameters copied.

deepcopy()[source]#

Return a deep copy of this model.

Returns:
resultGriddedPSFModel

A deep copy of this model.

evaluate(x, y, flux, x_0, y_0)[source]#

Calculate the ePSF model at the input coordinates for the given model parameters.

Parameters:
x, yfloat or ndarray

The x and y positions at which to evaluate the model.

fluxfloat

The flux scaling factor for the model.

x_0, y_0float

The (x, y) position of the model.

Returns:
evaluated_modelndarray

The evaluated model.

get_bounding_box(with_bbox=True)#

Return the bounding_box of a model if it exists or None otherwise.

Parameters:
with_bbox

The value of the with_bounding_box keyword argument when calling the model. Default is True for usage when looking up the model’s bounding_box without risk of error.

has_inverse()#

Returns True if the model has an analytic or user inverse defined.

input_shape(inputs)#

Get input shape for bounding_box evaluation.

output_units(**kwargs)#

Return a dictionary of output units for this model given a dictionary of fitting inputs and outputs.

The input and output Quantity objects should be given as keyword arguments.

Notes

This method is needed in order to be able to fit models with units in the parameters, since we need to temporarily strip away the units from the model during the fitting (which might be done by e.g. scipy functions).

This method will force extra model evaluations, which maybe computationally expensive. To avoid this, one can add a return_units property to the model, see return_units.

plot_grid(*, ax=None, vmax_scale=None, peak_norm=False, deltas=False, cmap='viridis', dividers=True, divider_color='darkgray', divider_ls='-', figsize=None)#

Plot the grid of ePSF models.

Parameters:
axmatplotlib.axes.Axes or None, optional

The matplotlib axes on which to plot. If None, then the current Axes instance is used.

vmax_scalefloat, optional

Scale factor to apply to the image stretch limits. This value is multiplied by the peak ePSF value to determine the plotting vmax. The defaults are 1.0 for plotting the ePSF data and 0.03 for plotting the ePSF difference data (deltas=True). If deltas=True, the vmin is set to -vmax. If deltas=False the vmin is set to vmax / 1e4.

peak_normbool, optional

Whether to normalize the ePSF data by the peak value. The default shows the ePSF flux per pixel.

deltasbool, optional

Set to True to show the differences between each ePSF and the average ePSF.

cmapstr or matplotlib.colors.Colormap, optional

The colormap to use. The default is ‘viridis’.

dividersbool, optional

Whether to show divider lines between the ePSFs.

divider_color, divider_lsstr, optional

Matplotlib color and linestyle options for the divider lines between ePSFs. These keywords have no effect unless show_dividers=True.

figsize(float, float), optional

The figure (width, height) in inches.

Returns:
figmatplotlib.figure.Figure

The matplotlib figure object. This will be the current figure if ax=None. Use fig.savefig() to save the figure to a file.

Notes

This method returns a figure object. If you are using this method in a script, you will need to call plt.show() to display the figure. If you are using this method in a Jupyter notebook, the figure will be displayed automatically.

When in a notebook, if you do not store the return value of this function, the figure will be displayed twice due to the REPL automatically displaying the return value of the last function call. Alternatively, you can append a semicolon to the end of the function call to suppress the display of the return value.

prepare_inputs(*inputs, model_set_axis=None, equivalencies=None, **kwargs)#

This method is used in __call__ to ensure that all the inputs to the model can be broadcast into compatible shapes (if one or both of them are input as arrays), particularly if there are more than one parameter sets. This also makes sure that (if applicable) the units of the input will be compatible with the evaluate method.

prepare_outputs(broadcasted_shapes, *outputs, **kwargs)#
classmethod rename(name=None, inputs=None, outputs=None)#

Return a copy of this model with a new name.

render(out=None, coords=None)#

Evaluate a model at fixed positions, respecting the bounding_box.

The key difference relative to evaluating the model directly is that this method is limited to a bounding box if the bounding_box attribute is set.

Parameters:
outnumpy.ndarray, optional

An array that the evaluated model will be added to. If this is not given (or given as None), a new array will be created.

coordsarray-like, optional

An array to be used to translate from the model’s input coordinates to the out array. It should have the property that self(coords) yields the same shape as out. If out is not specified, coords will be used to determine the shape of the returned array. If this is not provided (or None), the model will be evaluated on a grid determined by bounding_box.

Returns:
outnumpy.ndarray

The model added to out if out is not None, or else a new array from evaluating the model over coords. If out and coords are both None, the returned array is limited to the bounding_box limits. If bounding_box is None, arr or coords must be passed.

Raises:
ValueError

If coords are not given and the bounding_box of this model is not set.

Examples

Efficient Model Rendering with Bounding Boxes

set_slice_args(*args)#
strip_units_from_tree()#
sum_of_implicit_terms(*args, **kwargs)#

Evaluate the sum of any implicit model terms on some input variables. This includes any fixed terms used in evaluating a linear model that do not have corresponding parameters exposed to the user. The prototypical case is astropy.modeling.functional_models.Shift, which corresponds to a function y = a + bx, where b=1 is intrinsically fixed by the type of model, such that sum_of_implicit_terms(x) == x. This method is needed by linear fitters to correct the dependent variable for the implicit term(s) when solving for the remaining terms (ie. a = y - bx).

with_units_from_data(**kwargs)#

Return an instance of the model which has units for which the parameter values are compatible with the data units specified.

The input and output Quantity objects should be given as keyword arguments.

Notes

This method is needed in order to be able to fit models with units in the parameters, since we need to temporarily strip away the units from the model during the fitting (which might be done by e.g. scipy functions).

The units that the parameters will gain are not necessarily the units of the input data, but are derived from them. Model subclasses that want fitting to work in the presence of quantities need to define a _parameter_units_for_data_units method that takes the input and output units (as two dictionaries) and returns a dictionary giving the target units for each parameter.

without_units_for_data(**kwargs)#

Return an instance of the model for which the parameter values have been converted to the right units for the data, then the units have been stripped away.

The input and output Quantity objects should be given as keyword arguments.

Notes

This method is needed in order to be able to fit models with units in the parameters, since we need to temporarily strip away the units from the model during the fitting (which might be done by e.g. scipy functions).

The units that the parameters should be converted to are not necessarily the units of the input data, but are derived from them. Model subclasses that want fitting to work in the presence of quantities need to define a _parameter_units_for_data_units method that takes the input and output units (as two dictionaries) and returns a dictionary giving the target units for each parameter.