CompensableImage

class lsst.ts.wep.cwfs.CompensableImage(centroidFindType=CentroidFindType.RandomWalk)

Bases: object

Instantiate the class of CompensableImage.

Parameters:
centroidFindTypeenum ‘CentroidFindType’, optional

Algorithm to find the centroid of donut. (the default is CentroidFindType.RandomWalk.)

Methods Summary

centerOnProjection(img, template[, window])

Center the image to the template's center.

compensate(inst, algo, zcCol, model)

Calculate the image compensated from the affection of wavefront.

createBlendedCoadd(maskArray, blendPadding)

Cut out regions of the input mask where blended donuts overlap with the original donut.

getDefocalType()

Get the defocal type.

getFieldXY()

Get the field x, y in degree.

getImg()

Get the image.

getImgInit()

Get the initial image before doing the compensation.

getImgObj()

Get the image object.

getImgSizeInPix()

Get the image size in pixel.

getNonPaddedMask()

Get the non-padded mask corresponding to aperture.

getOffAxisCoeff()

Get the coefficients to do the off-axis correction.

getPaddedMask()

Get the padded mask use at the offset planes.

imageCoCenter(inst[, fov, debugLevel])

Shift the weighting center of donut to the center of reference image with the correction of projection of fieldX and fieldY. Parameters ---------- inst : Instrument Instrument to use. fov : float, optional Field of view (FOV) of telescope. (the default is 3.5.) debugLevel : int, optional Show the information under the running. If the value is higher, the information shows more. It can be 0, 1, 2, or 3. (the default is 0.).

isCaustic()

The image is caustic or not.

makeBlendedMask(inst, model, boundaryT, ...)

Create a binary mask of a central source with the area overlapping a blended object cut out.

makeMask(inst, model, boundaryT, ...[, ...])

Get the binary mask which considers the obscuration and off-axis correction.

makeMaskList(inst, model)

Calculate the mask list based on the obscuration and optical model.

setImg(fieldXY, defocalType[, blendOffsets, ...])

Set the wavefront image.

setOffAxisCorr(inst, order)

Set the coefficients of off-axis correction for x, y-projection of intra- and extra-image.

updateImage(image)

Update the image of donut.

updateImgInit()

Update the backup of initial image.

Methods Documentation

centerOnProjection(img, template, window=20)

Center the image to the template’s center.

Parameters:
imgnumpy.array

Image to be centered with the template. The input image needs to be a n-by-n matrix.

templatenumpy.array

Template image to have the same dimension as the input image (‘img’). The center of template is the position of input image tries to align with.

windowint, optional

Size of window in pixel. Assume the difference of centers of input image and template is in this range (e.g. [-window/2, window/2] if 1D). (the default is 20.)

Returns:
numpy.array

Recentered image.

compensate(inst, algo, zcCol, model)

Calculate the image compensated from the affection of wavefront.

Parameters:
instInstrument

Instrument to use.

algoAlgorithm

Algorithm to solve the Poisson’s equation. It can by done by the fast Fourier transform or serial expansion.

zcColnumpy.ndarray

Coefficients of wavefront.

modelstr

Optical model. It can be “paraxial”, “onAxis”, or “offAxis”.

Raises:
RuntimeError

input:size zcCol in compensate needs to be a numTerms row column vector.

createBlendedCoadd(maskArray, blendPadding)

Cut out regions of the input mask where blended donuts overlap with the original donut.

Parameters:
maskArraynumpy.ndarray

Original input binary mask with just a single unblended donut.

blendPaddingint

Number of pixels to increase the radius and expand the footprint of the blended source.

Returns:
numpy.ndarray [int]

New maskArray modified to remove the footprint of any overlapping donuts from the masked area of the original donut.

getDefocalType()

Get the defocal type.

Returns:
enum ‘DefocalType’

Defocal type.

getFieldXY()

Get the field x, y in degree.

Returns:
float

Field x in degree.

float

Field y in degree.

getImg()

Get the image.

Returns:
numpy.ndarray

Image.

getImgInit()

Get the initial image before doing the compensation.

Returns:
numpy.ndarray

Initial image before doing the compensation.

getImgObj()

Get the image object.

Returns:
Image

Image object.

getImgSizeInPix()

Get the image size in pixel.

Returns:
int

Image size in pixel.

getNonPaddedMask()

Get the non-padded mask corresponding to aperture.

Returns:
numpy.ndarray[int]

Non-padded mask

getOffAxisCoeff()

Get the coefficients to do the off-axis correction.

Returns:
numpy.ndarray

Coefficients to do the off-axis correction.

float

Defocused offset in files of off-axis correction.

getPaddedMask()

Get the padded mask use at the offset planes.

Returns:
numpy.ndarray[int]

Padded mask.

imageCoCenter(inst, fov=3.5, debugLevel=0)

Shift the weighting center of donut to the center of reference image with the correction of projection of fieldX and fieldY. Parameters ———- inst : Instrument

Instrument to use.

fovfloat, optional

Field of view (FOV) of telescope. (the default is 3.5.)

debugLevelint, optional

Show the information under the running. If the value is higher, the information shows more. It can be 0, 1, 2, or 3. (the default is 0.)

isCaustic()

The image is caustic or not.

The image might be caustic from the over-compensation.

Returns:
bool

True if the image is caustic.

makeBlendedMask(inst, model, boundaryT, maskScalingFactorLocal, blendPadding=8, compensated=True)

Create a binary mask of a central source with the area overlapping a blended object cut out. Thus this mask will cover only the unblended footprint of the central donut. Modifies self.mask_pupil and self.mask_comp (see makeMask for more details).

Parameters:
instInstrument

Instrument to use.

modelstr

Optical model. It can be “paraxial”, “onAxis”, or “offAxis”.

boundaryTint

Extended boundary in pixel. It defines how far the computation mask extends beyond the pupil mask. And, in fft, it is also the width of Neuman boundary where the derivative of the wavefront is set to zero.

maskScalingFactorLocalfloat

Mask scaling factor (for fast beam) for local correction.

blendPaddingint, optional

Number of pixels to increase the radius and expand the footprint of the blended source. (the default is 8.)

compensated: bool, optional

Whether to calculate masks for the compensated or original image. If True, mask will be in orientation of the compensated image. If False, mask will be in orientation of the original image. Note this is only relevant for extrafocal images. (the default is True.)

makeMask(inst, model, boundaryT, maskScalingFactorLocal, compensated=True)

Get the binary mask which considers the obscuration and off-axis correction.

There will be two mask parameters to be calculated: mask_comp: computation mask, i.e. padded mask,

for use at the offset planes

mask_pupil: pupil mask, i.e. non-padded mask,

corresponding to aperture

Parameters:
instInstrument

Instrument to use.

modelstr

Optical model. It can be “paraxial”, “onAxis”, or “offAxis”.

boundaryTint

Extended boundary in pixel. It defines how far the computation mask extends beyond the pupil mask. And, in fft, it is also the width of Neuman boundary where the derivative of the wavefront is set to zero.

maskScalingFactorLocalfloat

Mask scaling factor (for fast beam) for local correction.

compensated: bool, optional

Whether to calculate masks for the compensated or original image. If True, mask will be in orientation of the compensated image. If False, mask will be in orientation of the original image. Note this is only relevant for extrafocal images. (the default is True.)

makeMaskList(inst, model)

Calculate the mask list based on the obscuration and optical model.

Parameters:
instInstrument

Instrument to use.

modelstr

Optical model. It can be “paraxial”, “onAxis”, or “offAxis”.

Returns:
numpy.ndarray

The list of mask.

setImg(fieldXY, defocalType, blendOffsets=None, image=None, imageFile=None)

Set the wavefront image.

Parameters:
fieldXYtuple or list

Position of donut on the focal plane in degree (field x, field y).

defocalTypeenum ‘DefocalType’

Defocal type of image.

blendOffsetslist or None, optional

Positions of blended donuts relative to location of center donut. Enter as [xCoordList, yCoordList]. Length of xCoordList and yCoordList must be same length. (the default is None).

imagenumpy.ndarray, optional

Array of image. (the default is None.)

imageFilestr, optional

Path of image file. (the default is None.)

Raises:
ValueError

Input to blendOffsets must have same number of x-coordinates and y-coordinates.

setOffAxisCorr(inst, order)

Set the coefficients of off-axis correction for x, y-projection of intra- and extra-image.

This is for the mapping of coordinate from the telescope aperture to defocal image plane.

Parameters:
instInstrument

Instrument to use.

orderint

Up to order-th of off-axis correction.

updateImage(image)

Update the image of donut.

Parameters:
imagenumpy.ndarray

Donut image.

updateImgInit()

Update the backup of initial image.

This will be used in the outer loop iteration, which always uses the initial image (image0) before each iteration starts.