scarlet.interpolation

scarlet.interpolation.apply_2D_trapezoid_rule(y, x, f, dNy, dNx=None, dy=None, dx=None)[source]

Use the trapezoid rule to integrate over a subsampled function 2D implementation of the trapezoid rule. See apply_trapezoid_rule for a description, with the difference that f is a function f(y,x), where we note the c ++`(y,x)` ordering.

scarlet.interpolation.bilinear(dx)[source]

Bilinear interpolation kernel

Interpolate between neighboring pixels to shift by a fractional amount.

Parameters
dx: float

Fractional amount that the kernel will be shifted

Returns
result: array

2x2 linear kernel to use nearest neighbor interpolation

window: array

The pixel values for the window containing the kernel

scarlet.interpolation.catmull_rom(dx)[source]

Cubic spline with a=0.5, b=0

See cubic_spline for details.

scarlet.interpolation.common_projections(img1, img2)[source]

Project two images to a common frame

It is assumed that the two images have the same center. This is mainly used for FFT convolutions of source components, where the convolution kernel is a different size than the morphology.

Parameters
img1: array

1st 2D image to project

img2: array

2nd 2D image to project

Returns
img1: array

Projection of 1st image

img2: array

Projection of 2nd image

scarlet.interpolation.cubic_spline(dx, a=1, b=0)[source]

Generate a cubix spline centered on dx.

Parameters
dx: float

Fractional amount that the kernel will be shifted

a: float

Cubic spline sharpness paremeter

b: float

Cubic spline shape parameter

Returns
result: array

Cubic Spline kernel in a window from floor(dx)-1 to floor(dx) + 3

window: array

The pixel values for the window containing the kernel

scarlet.interpolation.get_affine(wcs)[source]
scarlet.interpolation.get_angles(frame_wcs, model_wcs)[source]

Computes the angles between two WCS Parameters ———-

frame_wcs: WCS

WCS of the observation’s frame

model_WCS:

WCS of the model frame.

scarlet.interpolation.get_common_padding(img1, img2, padding=None)[source]

Project two images to a common frame

It is assumed that the two images have the same center. This is mainly used for FFT convolutions of source components, where the convolution kernel is a different size than the morphology.

Parameters
img1: array

1st 2D or 3D image to project

img2: array

2nd 2D or 3D image to project

Returns
img1: array

Projection of 1st image

img2: array

Projection of 2nd image

scarlet.interpolation.get_filter_bounds(coords)[source]

Get the slices in x and y to apply a filter

Parameters
coords: array

The coordinates of the filter, defined by get_filter_coords.

Returns
y_start, y_end, x_start, x_end: int

The start and end of each slice that is passed to apply_filter.

scarlet.interpolation.get_filter_coords(filter_values, center=None)[source]

Create filter coordinate grid needed for the apply filter function

Parameters
filter_values: array

The 2D array of the filter to apply.

center: tuple

The center (y,x) of the filter. If center is None then filter_values must have an odd number of rows and columns and the center will be set to the center of filter_values.

Returns
coords: array

The coordinates of the pixels in filter_values, where the coordinates of the center pixel are (0,0).

scarlet.interpolation.get_pixel_size(model_affine)[source]

Extracts the pixel size from a wcs

scarlet.interpolation.get_projection_slices(image, shape, yx0=None)[source]

Get slices needed to project an image

This method returns the bounding boxes needed to project image into a larger image with shape. The results can be used with projection[bb] = image[ibb].

Parameters
image: array

2D input image

shape: tuple

Shape of the new image.

yx0: tuple

Location of the lower left corner of the image in the projection. If yx0 is None then the image is centered in the projection.

Returns
bb: tuple

(yslice, xslice) of the projected image to place image.

ibb: tuple

(iyslice, ixslice) of image to insert into the projection.

bounds: tuple

(bottom, top, left, right) locations of the corners of image in the projection. While this isn’t needed for slicing it can be useful for calculating information about the image before projection.

scarlet.interpolation.get_psf_size(psf)[source]

Measures the size of a psf by computing the size of the area in 3 sigma around the center.

This is an approximate method to estimate the size of the psf for setting the size of the frame, which does not require a precise measurement.

Parameters
PSF: `scarlet.PSF` object

PSF for whic to compute the size

Returns
——-
sigma3: float

radius of the area inside 3 sigma around the center in pixels

scarlet.interpolation.get_separable_kernel(dy, dx, kernel=<function lanczos>, **kwargs)[source]

Create a 2D kernel from a 1D kernel separable in x and y

Parameters
dy: float

amount to shift image in x

dx: float

amount to shift image in y

kernel: function

1D kernel that is separable in x and y

kwargs: dict

Keyword arguments for the kernel

Returns
kernel: Tensor

2D separable kernel

x_window: Tensor

The pixel values for the window containing the kernel in the x-direction

y_window: Tensor

The pixel values for the window containing the kernel in the y-direction

scarlet.interpolation.interpolate_observation(observation, frame, wave_filter=False)[source]

Interpolates the images in an observation on the grid described by a frame and returns the interpolated images.

Returns
interp: numpy.ndarray

array containing the interpolated images from observation.

scarlet.interpolation.lanczos(dx, a=3)[source]

Lanczos kernel

Parameters
dx: float

amount to shift image

a: int

Lanczos window size parameter

Returns
result: array-like

1D Lanczos kernel

scarlet.interpolation.mitchel_netravali(dx)[source]

Cubic spline with a=1/3, b=1/3

See cubic_spline for details.

scarlet.interpolation.mk_shifter(shape, real=False)[source]

Performs shifts in the Fourier domain on Fourier objects

scarlet.interpolation.project_image(image, shape, yx0=None)[source]

Project an image centered in a larger image

The projection pads the image with zeros if necessary or trims the edges if img is larger than shape in a given dimension.

Parameters
image: array

2D input image

shape: tuple

Shape of the new image.

yx0: tuple

Location of the lower left corner of the image in the projection. If yx0 is None then the image is centered in the projection.

Returns
result: array

The result of projecting image.

scarlet.interpolation.quintic_spline(dx, dtype=<class 'numpy.float64'>)[source]
scarlet.interpolation.sinc2D(y, x)[source]

2-D sinc function based on the product of 2 1-D sincs

Parameters
x, y: arrays

Coordinates where to evaluate the 2-D sinc

Returns
——-
result: array

2-D sinc evaluated in x and y

scarlet.interpolation.sinc_interp(images, coord_hr, coord_lr, angle=None, padding=3)[source]
Parameters
image: array

image whose pixels are at positions coord_lr

coord_hr: array (2xN)

Coordinates of the high resolution grid

coord_lr: array (2xM)

Coordinates of the low resolution grid

angle: float

rotation angle between coordinate sets coord_hr and coord_lr

padding: int

value of zero padding for fft

Returns
——-

result: interpolated samples at positions coord_hr

scarlet.interpolation.sinc_interp_inplace(image, h_image, h_target, angle, pad_shape=None)[source]

In place interpolation of a cube of images

Performs interpolation from a grid defined by the grid of image to a grid spanning the same physical area scaled by a factor h and rotated by angle radians. The center for the rotation is the central pixel of the image. This procedure is advised for odd-sized images only.

Parameters
image: `ndarray`

Cube of images with shape BxNyxNx with B the number of bands and NyxNx, the number of pixels

h_image: `float`

Phisical scale of a pixel in image

h_target: `float`

Physical scale of the target pixel to which to interpolate

angle: float

angle between the grid of image and the target grid where to interpolate

Returns
——-
interp_image: `ndarray`

padded interpolated image

scarlet.interpolation.subsample_function(y, x, f, dNy, dNx=None, dy=None, dx=None)[source]

Subsample a function Given the expected pixel grid of a function, subsample that function at a grid subdivided in x by dNx and y by dNy.