psfmodel Module

Point-spread function (PSF) models and fitting for image data.

This package defines an abstract PSF API and concrete implementations such as GaussianPSF for analytic Gaussians, pixel integration, background handling, and astrometric fitting.

The public surface re-exports PSF and GaussianPSF (see __all__). Additional modules (for example instrument-specific PSFs) are imported from their submodules when needed. A logging.NullHandler is attached to the package logger so library logging is opt-in for applications.

class PSF(*, logger: Logger | None = None, detailed_logging: bool = False, **kwargs: Any)[source]

Bases: ABC

Abstract base for 2-D point-spread models used in fitting and rendering.

Subclass PSF to provide a concrete model. The base class supplies shared utilities (for example find_position(), background helpers, and motion smear via _eval_rect_smeared()); evaluation itself is defined by subclasses through the abstract API below.

Abstract methods (must implement)

  • eval_point() – Evaluate the continuous PSF at fractional (y, x) (scalar or broadcast arrays). Signature includes keyword-only scale and base. Returns a float or a numpy.ndarray of floats matching the broadcast shape of coord. Origin (0, 0) is the PSF center; coordinates may be negative. Subclasses may add keyword-only parameters.

  • eval_rect() – Build a rectangular, pixel-integrated patch. Signature: rect_size, offset, then keyword-only movement, movement_granularity, scale, base, and subclass-specific **kwargs. Must return numpy.ndarray with dtype float64 and shape (height, width) matching rect_size as (size_y, size_x). Should validate inputs (for example odd rect_size) and raise ValueError for invalid arguments with a clear message.

  • _eval_rect() – Internal hook for the same patch without the checks in eval_rect(); same core keyword-only scale and base. Must return a float64 array of shape (height, width). Subclasses often add keyword-only model parameters. Callers must pass consistent arguments; implementations may omit validation.

Protected attributes

  • _loggerlogging.Logger for this instance (set in __init__()). Subclasses log warnings and diagnostics through it.

  • _additional_paramslist (initialized empty), each entry a (lower_bound, upper_bound, name) tuple of two floats and a str keyword name. Used by find_position() / _find_position() to append extra optimized parameters (bounds and eval_rect() keyword). Subclasses append one tuple per fittable quantity in construction order; leave the list empty if there are no extra parameters (for example fixed-width models).

Errors and return conventions

Public evaluators should reject bad arguments with ValueError (or TypeError for wrong types) where feasible. Some higher-level routines (notably find_position()) signal failure by returning None instead of raising. Numeric outputs are real floating point; scale multiplies the model amplitude and base adds a constant offset in the same units as the evaluated PSF values unless a subclass documents physical units.

__init__(*, logger: Logger | None = None, detailed_logging: bool = False, **kwargs: Any) None[source]

Create a PSF object. Only called by subclasses.

Parameters:
  • logger – Logger for diagnostic messages from this instance. If omitted, uses logging.getLogger() with this module’s __name__.

  • detailed_logging – If True, emit INFO and DEBUG messages during fitting (for example from find_position()). Optimizer failure still logs at WARNING when this is False.

  • **kwargs – Reserved for subclass forward-compatibility; unknown names are ignored.

abstractmethod eval_point(coord: tuple[float | ndarray[tuple[Any, ...], dtype[floating]], float | ndarray[tuple[Any, ...], dtype[floating]]] | ndarray[tuple[Any, ...], dtype[floating]], *, scale: float = 1.0, base: float = 0.0) float | ndarray[tuple[Any, ...], dtype[floating]][source]

Evaluate the PSF at a single, fractional, point.

(0, 0) is the center of the PSF and x and y may be negative.

Parameters:
  • coord – The coordinate (y, x) at which to evaluate the PSF.

  • scale – A scale factor to apply to the resulting PSF.

  • base – A scalar added to the resulting PSF.

Other parameters may be available for specific subclasses.

Returns:

The PSF value at the given coordinate.

abstractmethod eval_rect(rect_size: list[int] | tuple[int, int], offset: list[float] | tuple[float, float] = (0.5, 0.5), *, movement: tuple[float, float] | None = None, movement_granularity: float = 0.1, scale: float = 1.0, base: float = 0.0, **kwargs: Any) ndarray[tuple[Any, ...], dtype[float64]][source]

Create a rectangular pixelated PSF.

This is done by evaluating the PSF function from [-rect_size_y//2:rect_size_y//2] to [-rect_size_x//2:rect_size_x//2].

The returned array has the PSF offset from the center by (offset_y, offset_x). An offset of (0, 0) places the PSF in the upper left corner of the center pixel while an offset of (0.5, 0.5) places the PSF in the center of the center pixel. The angle is applied relative to this new origin, so as angle changes the center of the ellipse does not move.

Parameters:
  • rect_size – The size of the rectangle (rect_size_y, rect_size_x) of the returned PSF. Both dimensions must be odd.

  • offset – The amount (offset_y, offset_x) to offset the center of the PSF.

  • movement – The amount of motion blur in the (Y, X) direction. Must be a tuple of scalars. None means no movement.

  • movement_granularity – The number of pixels to step for each smear while doing motion blur. A smaller granularity means that the resulting PSF will be more precise but also take longer to compute.

  • scale – A scale factor to apply to the resulting PSF.

  • base – A scalar added to the resulting PSF.

Other inputs may be available for specific subclasses.

Returns:

The integral of the 2-D PSF over each full pixel in the rectangle.

static background_gradient_fit(image: ndarray[tuple[Any, ...], dtype[floating]], order: int = 2, ignore_center: int | tuple[int, int] | None = None, num_sigma: float | None = None, debug: bool = False, *, logger: Logger | None = None) tuple[ndarray[tuple[Any, ...], dtype[float64]] | None, ndarray[tuple[Any, ...], dtype[bool]] | None][source]

Return the polynomial fit to the pixels of an image.

Parameters:
  • image – 2D array to fit; must have odd shape in each dimension.

  • order – Order of the polynomial.

  • ignore_center – A scalar or tuple (ignore_y, ignore_x) giving the number of pixels on either side of the center to ignore while fitting. 0 means ignore the center pixel. None means don’t ignore anything.

  • num_sigma – Outlier rejection uses the fit residual image - gradient: unmasked pixels with absolute residual at least num_sigma times the standard deviation of that residual (mask-aware) are masked and the fit is repeated until convergence or until sigma falls below the numerical noise floor. None disables this. Non-positive values disable masking after the initial least-squares fit.

  • debug – Set to debug bad pixel removal.

  • logger – Logger for debug messages; defaults to this module’s logger.

Returns:

A tuple of the background coefficient array and the mask of ignored pixels.

static background_gradient(rect_size: tuple[int, int], bkgnd_params: ArrayLike) ndarray[tuple[Any, ...], dtype[float64]][source]

Create a background gradient.

Parameters:
  • rect_size(size_y, size_x), the shape of the output grid (height, width) in pixels; must match the image shape used when the coefficients were fit.

  • bkgnd_params – Coefficients of the background polynomial (1-D array-like). The polynomial order is inferred from the number of elements.

Returns:

A numpy.ndarray of dtype float64 with shape rect_size (i.e. (size_y, size_x)): the evaluated 2-D background polynomial at each pixel center of the grid.

find_position(image: ndarray[tuple[Any, ...], dtype[floating]], box_size: tuple[int, int], starting_point: tuple[float, float], *, search_limit: float | tuple[float, float] = (1.5, 1.5), bkgnd_degree: int | None = 2, bkgnd_ignore_center: tuple[int, int] = (2, 2), bkgnd_num_sigma: float | None = None, tolerance: float = 1e-06, num_sigma: float | None = None, max_bad_frac: float = 0.2, allow_nonzero_base: bool = False, scale_limit: float = 1000.0, use_angular_params: bool = True, compute_uncertainty: bool = True) None | tuple[float, float, dict[str, Any]][source]

Find the (y, x) coordinates that best fit a 2-D PSF to an image.

Parameters:
  • image – The image (2-D).

  • box_size – A tuple (box_y, box_x) indicating the size of the PSF to use. This governs both the size of the PSF created at each step as well as the size of the subimage looked at. Both box_y and box_x must be odd.

  • starting_point – A tuple (y, x) indicating the best guess for where the object can be found. Searching is limited to a region around this point controlled by search_limit.

  • search_limit – A scalar or tuple (y_limit, x_limit) specifying the maximum distance to search from starting_point. If a scalar, both x_limit and y_limit are the same.

  • bkgnd_degree – The degree (order) of the background gradient polynomial. None means no background gradient is fit.

  • bkgnd_ignore_center – A tuple (ny, nx) giving the number of pixels on each side of the center point to ignore when fitting the background. The ignored region is thus ny*2+1 by nx*2+1.

  • bkgnd_num_sigma – The number of sigma a pixel needs to be beyond the background gradient to be ignored. None means don’t ignore bad pixels while computing the background gradient.

  • tolerance – The tolerance (both X and Function) in the Powell optimization algorithm.

  • num_sigma – The number of sigma for a pixel to be considered bad during PSF fitting. None means don’t ignore bad pixels while fitting the PSF.

  • max_bad_frac – The maximum allowable number of pixels masked during PSF fitting. If more pixels than this fraction are masked, the position fit fails.

  • allow_nonzero_base – If True, allow the base of the PSF (constant bias) to vary. Otherwise the base of the PSF is always at zero and can only scale in the positive direction.

  • scale_limit – The maximum PSF scale allowed.

  • use_angular_params – Use angles to optimize parameter values.

  • compute_uncertainty – If True (default), compute 1-sigma parameter uncertainties via a finite-difference Jacobian after the fit. This requires one additional forward-model evaluation per free parameter, roughly doubling total cost for a typical 3-parameter fit. Set to False in hot loops (e.g. batch characterization runs) when the x_err, y_err, scale_err, and base_err metadata entries are not needed; they will be NaN when skipped.

Returns:

None if no fit found.

Otherwise returns pos_y, pos_x, metadata. Metadata is a dictionary containing:

'x'                    Full-image X coordinate of fitted position
                       (same as pos_x).
'x_err'                1-sigma uncertainty in X (pixels).
'y'                    Full-image Y coordinate of fitted position
                       (same as pos_y).
'y_err'                1-sigma uncertainty in Y (pixels).
'_local_x'             Subimage-relative X offset (offset from the
                       center of the cropped subimage).
'_local_y'             Subimage-relative Y offset.
'scale'                The best fit PSF scale.
'scale_err'            1-sigma uncertainty in PSF scale.
'base'                 The best fit PSF base.
'base_err'             1-sigma uncertainty in PSF base; 0.0
                       when allow_nonzero_base is False.
'residual_rss'         Sum of squared residuals over unmasked
                       pixels.
'reduced_chi2'         residual_rss divided by degrees of
                       freedom (n_valid - n_params); near 1.0
                       for a noise-limited fit.
'noise_rms'            Per-pixel noise estimate from residuals:
                       sqrt(rss / n_valid).
'peak_snr'             Amplitude signal-to-noise ratio:
                       scale / noise_rms.
'subimg'               The box_size area of the original image
                       surrounding starting_point masked as
                       necessary using the num_sigma threshold.
'bkgnd_params'         The tuple of parameters defining the
                       background gradient.
'bkgnd_mask'           The mask used during background gradient
                       fitting.
'gradient'             The box_size background gradient.
'subimg-gradient'      The subimg with the background gradient
                       subtracted.
'psf'                  The PSF model from eval_rect with the
                       fitted scale and base (same array as
                       scaled_psf).
'scaled_psf'           Same as psf; model to compare to
                       subimg-gradient during outlier rejection.

In addition, metadata includes two entries for each “additional parameter” used during optimization: one for the value and one for the 1-sigma uncertainty ('param' and 'param_err').

Raises:
  • ValueError – If box_size is not a tuple of odd positive integers.

  • TypeError – If num_sigma is not a number or None.

  • ValueError – If num_sigma is not greater than 0.

  • ValueError – If the starting point is too close to the edge of the image.

  • ValueError – If the subimage has too many pixels masked.

class GaussianPSF(*, sigma: float | tuple[float | None, float | None] | None = None, mean: float | tuple[float, float] = 0.0, angle: float | None = 0.0, sigma_x_range: tuple[float, float] = (0.01, 10.0), sigma_y_range: tuple[float, float] = (0.01, 10.0), angle_subsample: int = 13, logger: Logger | None = None, detailed_logging: bool = False)[source]

Bases: PSF

A 2-D Gaussian symmetric PSF.

The PSF can have different standard deviations in the X and Y directions. The standard deviations for X and Y can be locked up front when the GaussianPSF object is created or left to float so that future calls may specify them directly.

Because these are so fast and easy to compute, we don’t cache any results.

__init__(*, sigma: float | tuple[float | None, float | None] | None = None, mean: float | tuple[float, float] = 0.0, angle: float | None = 0.0, sigma_x_range: tuple[float, float] = (0.01, 10.0), sigma_y_range: tuple[float, float] = (0.01, 10.0), angle_subsample: int = 13, logger: Logger | None = None, detailed_logging: bool = False) None[source]

Create a GaussianPSF object describing a 2-D Gaussian PSF.

Parameters:
  • sigma – The standard deviation of the Gaussian. May be a scalar, in which case the value applies to both X and Y, or a tuple (sigma_y, sigma_x) one of which may be None. None for sigma or for one of sigma_x/y means that the sigma will be supplied later.

  • mean – The mean of the Gaussian. May be a scalar in which case the value applies to both X and Y, or a tuple (mean_y, mean_x).

  • angle – Rotation angle from 0 to pi (0 is +X). Default 0.0 fixes an unrotated Gaussian and keeps pixel integration on the fast axis-aligned path. Pass None to include angle in _additional_params so it can be fitted in find_position(); per-call overrides are still accepted by eval_rect().

  • sigma_x_range – The valid range for sigma_x if it is not specified otherwise. This is used during PSF fitting to let sigma_x float to its optimal value.

  • sigma_y_range – The valid range for sigma_y if it is not specified otherwise. This is used during PSF fitting to let sigma_y float to its optimal value.

  • angle_subsample – The amount of subsampling to do in X and Y when computing a 2-D Gaussian pixel with a non-zero angle.

  • logger – Optional logger passed to PSF; see PSF.__init__().

  • detailed_logging – Passed to PSF; see PSF.__init__().

property sigma_y: float | None

Standard deviation of the Gaussian along the y-axis (pixels), if fixed.

Returns:

The locked sigma_y from construction, or None when it is left free for fitting or per-call overrides.

property sigma_x: float | None

Standard deviation of the Gaussian along the x-axis (pixels), if fixed.

Returns:

The locked sigma_x from construction, or None when it is left free for fitting or per-call overrides.

property mean_y: float

Center of the Gaussian along the y-axis in pixel coordinates.

Returns:

The mean y used when evaluating the model (scalar).

property mean_x: float

Center of the Gaussian along the x-axis in pixel coordinates.

Returns:

The mean x used when evaluating the model (scalar).

static gaussian_1d(x: float | ndarray[tuple[Any, ...], dtype[floating]], *, sigma: float = 1.0, mean: float = 0.0, scale: float = 1.0, base: float = 0.0) float | ndarray[tuple[Any, ...], dtype[floating]][source]

Return a 1-D Gaussian.

This simply returns the value of the Gaussian at a point or series of points. The Gaussian is normalized so the area under the curve is “scale”.

Parameters:
  • x – A scalar or array at which to evaluate the Gaussian function.

  • sigma – The standard deviation of the Gaussian.

  • mean – The mean of the Gaussian.

  • scale – The scale of the Gaussian; the area under the complete curve (excluding the base).

  • base – The base of the Gaussian; a scalar added to the curve.

Returns: The value of the Gaussian function at the point(s) defined by x.

static gaussian_2d(y: float | ndarray[tuple[Any, ...], dtype[floating]], x: float | ndarray[tuple[Any, ...], dtype[floating]], *, sigma_y: float = 1.0, sigma_x: float = 1.0, mean_y: float = 0.0, mean_x: float = 0.0, scale: float = 1.0, base: float = 0.0, angle: float = 0.0) float | ndarray[tuple[Any, ...], dtype[floating]][source]

Return a 2-D Gaussian using angle (angle 0-pi, 0 at 3 o’clock, CW).

This simply returns the value of the 2-D Gaussian at a series of points. The Gaussian is normalized so the area under the curve is scale.

Parameters:
  • y – A scalar or array at which to evaluate the 2-D Gaussian function.

  • x – A scalar or array at which to evaluate the 2-D Gaussian function.

  • sigma_y – The standard deviation of the Gaussian in the Y direction.

  • sigma_x – The standard deviation of the Gaussian in the X direction.

  • mean_y – The mean of the Gaussian in the Y dimension.

  • mean_x – The mean of the Gaussian in the X dimension.

  • scale – The scale of the Gaussian; the area under the complete curve (excluding the base).

  • base – The base of the Gaussian; a scalar added to the curve.

  • angle – The angle of the ellipse. Angle ranges from 0 to pi, with 0 being “3 o’clock” (+X) assuming that (0, 0) is in the top left corner. mean_y/x are specified in the rotated coordinate system.

Returns:

The value of the 2-D Gaussian at the points defined by x and y.

static gaussian_integral_1d(x_min: float | ndarray[tuple[Any, ...], dtype[floating]], x_max: float | ndarray[tuple[Any, ...], dtype[floating]], *, sigma: float = 1.0, mean: float = 0.0, scale: float = 1.0, base: float = 0.0) float | ndarray[tuple[Any, ...], dtype[floating]][source]

Return the integral of a Gaussian.

The integral is over the limits [xmin, xmax].

Values are generated via the error function, where the integral from -inf to x is proportional to

1 + erf((x - mean) / (sqrt(2) * sigma))

This function works for both scalar and array values of xmin and xmax.

Parameters:
  • x_min – The lower bound of the integral.

  • x_max – The upper bound of the integral.

  • sigma – The standard deviation of the Gaussian.

  • mean – The mean of the Gaussian.

  • scale – The scale of the Gaussian; the area under the complete curve (excluding the base).

  • base – The base of the Gaussian; a scalar added to the curve.

Returns:

The integral of the Gaussian between x_min and x_max.

Raises:

ValueError – If sigma is not positive.

static gaussian_integral_2d(y_min: float | ndarray[tuple[Any, ...], dtype[floating]], y_max: float | ndarray[tuple[Any, ...], dtype[floating]], x_min: float | ndarray[tuple[Any, ...], dtype[floating]], x_max: float | ndarray[tuple[Any, ...], dtype[floating]], *, sigma_y: float = 1.0, sigma_x: float = 1.0, mean_y: float = 0.0, mean_x: float = 0.0, scale: float = 1.0, base: float = 0.0, angle: float = 0.0, angle_subsample: int = 13) float | ndarray[tuple[Any, ...], dtype[floating]][source]

Return the double integral of a 2-D Gaussian.

The integral is over the limits [y_min, y_max] and [x_min, x_max].

For angle == 0 the integral is computed analytically using the error function (fast, exact). For a non-zero angle, the integral is approximated by uniform subsampling within each pixel using angle_subsample points per axis. When the bounds are arrays (e.g. a full image patch), all pixels are evaluated simultaneously in a single vectorised NumPy call on an (N, S, S) grid, where N is the number of pixels and S is angle_subsample; this avoids a per-pixel Python loop and is typically 10-20x faster than an equivalent element-wise approach.

Parameters:
  • y_min – The lower bound of the integral in the Y dimension.

  • y_max – The upper bound of the integral in the Y dimension.

  • x_min – The lower bound of the integral in the X dimension.

  • x_max – The upper bound of the integral in the X dimension.

  • sigma_y – The standard deviation of the Gaussian in the Y dimension.

  • sigma_x – The standard deviation of the Gaussian in the X dimension.

  • mean_y – The mean of the Gaussian in the Y dimension.

  • mean_x – The mean of the Gaussian in the X dimension.

  • scale – The scale of the Gaussian; the area under the complete curve (excluding the base).

  • base – The base of the Gaussian; a scalar added to the curve.

  • angle – The angle of the ellipse. Angle ranges from 0 to pi, with 0 being “3 o’clock” (+X) assuming that (0, 0) is in the top left corner. Thus specifying 0 means that the X and Y axes of the ellipse align with the standard X and Y Cartesian axes. mean_y/x are specified in the rotated coordinate system, while sigma_y/x are specified in the unrotated Cartesian coordinate system.

  • angle_subsample – The amount of subsampling to do in X and Y when computing a 2-D Gaussian pixel with a non-zero angle.

Returns:

The integral of the 2-D Gaussian between y_min and y_max, and x_min and x_max.

eval_point(coord: tuple[float | ndarray[tuple[Any, ...], dtype[floating]], float | ndarray[tuple[Any, ...], dtype[floating]]] | ndarray[tuple[Any, ...], dtype[floating]], *, sigma: float | None = None, scale: float = 1.0, base: float = 0.0, sigma_y: float | None = None, sigma_x: float | None = None, angle: float | None = None) float | ndarray[tuple[Any, ...], dtype[floating]][source]

Evaluate the 2-D Gaussian PSF at a single, fractional, point.

(0, 0) is the center of the PSF and x and y may be negative.

Parameters:
  • coord – The coordinate (y, x) at which to evaluate the PSF.

  • scale – A scale factor to apply to the resulting PSF.

  • base – A scalar added to the resulting PSF.

  • sigma – Standard deviations: a scalar (both axes) or a (sigma_y, sigma_x) tuple. Must not duplicate values fixed at construction; None if both are already set on the instance.

  • sigma_y – An alternative way to specify sigma_y. Used primarily for letting sigma_y float during PSF fitting.

  • sigma_x – An alternative way to specify sigma_x. Used primarily for letting sigma_x float during PSF fitting.

  • angle – An alternative way to specify angle. Used primarily for letting angle float during PSF fitting.

Returns:

The value of the 2-D Gaussian at the point(s) specified by coord.

eval_pixel(coord: tuple[int | ndarray[tuple[Any, ...], dtype[int64]], int | ndarray[tuple[Any, ...], dtype[int64]]] | ndarray[tuple[Any, ...], dtype[floating]], offset: tuple[float | ndarray[tuple[Any, ...], dtype[floating]], float | ndarray[tuple[Any, ...], dtype[floating]]] | ndarray[tuple[Any, ...], dtype[floating]] = (0.5, 0.5), *, scale: float = 1.0, base: float = 0.0, sigma: tuple[float, float] | None = None, sigma_y: float | None = None, sigma_x: float | None = None, angle: float | None = None) float | ndarray[tuple[Any, ...], dtype[floating]][source]

Evaluate the Gaussian PSF integrated over an entire integer pixel.

The returned array has the PSF offset from the center by (offset_y, offset_x). An offset of (0, 0) places the PSF in the upper left corner of the center pixel while an offset of (0.5, 0.5) places the PSF in the center of the center pixel. The angle is applied relative to this new origin, so as angle changes the center of the ellipse does not move.

This essentially performs a 2-D integration of the PSF over the intervals [y-offset_y-0.5, y-offset_y+0.5] and [x-offset_x-0.5, x-offset_x+0.5].

Parameters:
  • coord – The integer coordinate(s) (y, x) at which to evaluate the PSF.

  • offset – The amount (offset_y, offset_x) to offset the center of the PSF.

  • scale – A scale factor to apply to the resulting PSF.

  • base – A scalar added to the resulting PSF.

  • sigma – Standard deviations: a scalar (both axes) or a (sigma_y, sigma_x) tuple. Must not duplicate values fixed at construction; None if both are already set on the instance.

  • sigma_y – An alternative way to specify sigma_y. Used primarily for letting sigma_y float during PSF fitting.

  • sigma_x – An alternative way to specify sigma_x. Used primarily for letting sigma_x float during PSF fitting.

  • angle – An alternative way to specify angle. Used primarily for letting angle float during PSF fitting.

Returns:

The integral of the 2-D Gaussian over each full pixel.

eval_rect(rect_size: list[int] | tuple[int, int], offset: list[float] | tuple[float, float] = (0.5, 0.5), *, movement: tuple[float, float] | None = None, movement_granularity: float = 0.1, scale: float = 1.0, base: float = 0.0, sigma: tuple[float, float] | None = None, sigma_y: float | None = None, sigma_x: float | None = None, angle: float | None = None, **kwargs: Any) ndarray[tuple[Any, ...], dtype[float64]][source]

Create a rectangular pixelated Gaussian PSF.

This is done by evaluating the PSF function from [-rect_size_y//2:rect_size_y//2] to [-rect_size_x//2:rect_size_x//2].

The returned array has the PSF offset from the center by (offset_y, offset_x). An offset of (0, 0) places the PSF in the upper left corner of the center pixel while an offset of (0.5, 0.5) places the PSF in the center of the center pixel. The angle is applied relative to this new origin, so as angle changes the center of the ellipse does not move.

Parameters:
  • rect_size – The size of the rectangle (rect_size_y, rect_size_x) of the returned PSF. Both dimensions must be odd.

  • offset – The amount (offset_y, offset_x) to offset the center of the PSF.

  • movement – The amount of motion blur in the (Y, X) direction. Must be a tuple of scalars. None means no movement.

  • movement_granularity – The number of pixels to step for each smear while doing motion blur. A smaller granularity means that the resulting PSF will be more precise but also take longer to compute.

  • scale – A scale factor to apply to the resulting PSF.

  • base – A scalar added to the resulting PSF.

  • sigma – Standard deviations: a scalar (both axes) or a (sigma_y, sigma_x) tuple. Must not duplicate values fixed at construction; None if both are already set on the instance.

  • sigma_y – An alternative way to specify sigma_y. Used primarily for letting sigma_y float during PSF fitting.

  • sigma_x – An alternative way to specify sigma_x. Used primarily for letting sigma_x float during PSF fitting.

  • angle – An alternative way to specify angle. Used primarily for letting angle float during PSF fitting.

Returns:

The integral of the 2-D Gaussian over each full pixel in the rectangle.