autocti.Mask2D#

class autocti.Mask2D(mask: Union[ndarray, List], pixel_scales: Union[Tuple[float], Tuple[float, float], float], sub_size: int = 1, origin: Tuple[float, float] = (0.0, 0.0), invert: bool = False, *args, **kwargs)[source]#

Bases: Mask2D

A 2D mask, used for masking values which are associated with a a uniform rectangular grid of pixels.

When applied to 2D data with the same shape, values in the mask corresponding to False entries are unmasked and therefore used in subsequent calculations. .

The ``Mask2D`, has in-built functionality which:

  • Maps data structures between two data representations: slim` (all unmasked False values in a 1D ndarray) and native (all unmasked values in a 2D or 3D ndarray).

  • Has a Geometry2D object (defined by its (y,x) pixel scales, (y,x) origin and sub_size) which defines how coordinates are converted from pixel units to scaled units.

  • Associates Cartesian Grid2D objects of (y,x) coordinates with the data structure (e.g. a (y,x) grid of all unmasked pixels) via the DeriveGrid2D object.

  • This includes sub-grids, which perform calculations higher resolutions which are then binned up.

A detailed description of the 2D mask API is provided below.

SLIM DATA REPRESENTATION (sub-size=1)

Below is a visual illustration of a Mask2D, where a total of 10 pixels are unmasked (values are False):

x x x x x x x x x x
x x x x x x x x x x     This is an example ``Mask2D``, where:
x x x x x x x x x x
x x x x O O x x x x     x = `True` (Pixel is masked and excluded from the array)
x x x O O O O x x x     O = `False` (Pixel is not masked and included in the array)
x x x O O O O x x x
x x x x x x x x x x
x x x x x x x x x x
x x x x x x x x x x
x x x x x x x x x x

The mask pixel index’s are as follows (the positive / negative direction of the Grid2D objects associated with the mask are also shown on the y and x axes).

<--- -ve  x  +ve -->

 x x x x x x x x x x  ^
 x x x x x x x x x x  I
 x x x x x x x x x x  I
 x x x x 0 1 x x x x +ve
 x x x 2 3 4 5 x x x  y
 x x x 6 7 8 9 x x x -ve
 x x x x x x x x x x  I
 x x x x x x x x x x  I
 x x x x x x x x x x \/
 x x x x x x x x x x

The Mask2D’s slim data representation is an ndarray of shape [total_unmasked_pixels].

For the Mask2D above the slim representation therefore contains 10 entries and two examples of these entries are:

mask[3] = the 4th unmasked pixel's value.
mask[6] = the 7th unmasked pixel's value.

A Cartesian grid of (y,x) coordinates, corresponding to all slim values (e.g. unmasked pixels) is given by mask.derive_grid.masked.slim.

NATIVE DATA REPRESENTATION (sub_size=1)

Masked data represented as an an ndarray of shape [total_y_values, total_x_values], where all masked entries have values of 0.0.

For the following mask:

x x x x x x x x x x
x x x x x x x x x x     This is an example ``Mask2D``, where:
x x x x x x x x x x
x x x x O O x x x x     x = `True` (Pixel is masked and excluded from the array)
x x x O O O O x x x     O = `False` (Pixel is not masked and included in the array)
x x x O O O O x x x
x x x x x x x x x x
x x x x x x x x x x
x x x x x x x x x x
x x x x x x x x x x

The mask has the following indexes:

<--- -ve  x  +ve -->

 x x x x x x x x x x  ^
 x x x x x x x x x x  I
 x x x x x x x x x x  I
 x x x x 0 1 x x x x +ve
 x x x 2 3 4 5 x x x  y
 x x x 6 7 8 9 x x x -ve
 x x x x x x x x x x  I
 x x x x x x x x x x  I
 x x x x x x x x x x  \/
 x x x x x x x x x x

In the above array:

- mask[0,0] = True (it is masked)
- mask[0,0] = True (it is masked)
- mask[3,3] = True (it is masked)
- mask[3,3] = True (it is masked)
- mask[3,4] = False (not masked)
- mask[3,5] = False (not masked)
- mask[4,5] = False (not masked)

SLIM TO NATIVE MAPPING

The Mask2D has functionality which maps data between the slim and native data representations.

For the example mask above, the 1D ndarray given by mask.derive_indexes.slim_to_native is:

slim_to_native[0] = [3,4]
slim_to_native[1] = [3,5]
slim_to_native[2] = [4,3]
slim_to_native[3] = [4,4]
slim_to_native[4] = [4,5]
slim_to_native[5] = [4,6]
slim_to_native[6] = [5,3]
slim_to_native[7] = [5,4]
slim_to_native[8] = [5,5]
slim_to_native[9] = [5,6]

SUB GRIDDING

If the Mask2D sub_size is > 1, its slim and native data representations have entries corresponding to the values at the centre of every sub-pixel of each unmasked pixel.

The sub-array indexes are ordered such that pixels begin from the first (top-left) sub-pixel in the first unmasked pixel. Indexes then go over the sub-pixels in each unmasked pixel, for every unmasked pixel.

Therefore, the shapes of the sub-array are as follows:

  • slim representation: an ndarray of shape [total_unmasked_pixels*sub_size**2].

  • native representation: an ndarray of shape [total_y_values*sub_size, total_x_values*sub_size].

Below is a visual illustration of a sub array. Indexing of each sub-pixel goes from the top-left corner. In contrast to the array above, our illustration below restricts the mask to just 2 pixels, to keep the illustration brief.

x x x x x x x x x x
x x x x x x x x x x     This is an example ``Mask2D``, where:
x x x x x x x x x x
x x x x x x x x x x     x = `True` (Pixel is masked and excluded from lens)
x 0 0 x x x x x x x     O = `False` (Pixel is not masked and included in lens)
x x x x x x x x x x
x x x x x x x x x x
x x x x x x x x x x
x x x x x x x x x x
x x x x x x x x x x

If sub_size=2, each unmasked pixel has 4 (2x2) sub-pixel values. For the example above, pixels 0 and 1 each have 4 values which map to slim representation as follows:

Pixel 0 - (2x2):

       slim[0] = value of first sub-pixel in pixel 0.
 0 1   slim[1] = value of first sub-pixel in pixel 1.
 2 3   slim[2] = value of first sub-pixel in pixel 2.
       slim[3] = value of first sub-pixel in pixel 3.

Pixel 1 - (2x2):

       slim[4] = value of first sub-pixel in pixel 0.
 4 5   slim[5] = value of first sub-pixel in pixel 1.
 6 7   slim[6] = value of first sub-pixel in pixel 2.
       slim[7] = value of first sub-pixel in pixel 3.

For the native data representation we get the following mappings:

Pixel 0 - (2x2):

       native[8, 2] = value of first sub-pixel in pixel 0.
 0 1   native[8, 3] = value of first sub-pixel in pixel 1.
 2 3   native[9, 2] = value of first sub-pixel in pixel 2.
       native[9, 3] = value of first sub-pixel in pixel 3.

Pixel 1 - (2x2):

       native[10, 4] = value of first sub-pixel in pixel 0.
 4 5   native[10, 5] = value of first sub-pixel in pixel 1.
 6 7   native[11, 4] = value of first sub-pixel in pixel 2.
       native[11, 5] = value of first sub-pixel in pixel 3.

Other entries (all masked sub-pixels are zero):

       native[0, 0] = 0.0 (it is masked, thus zero)
       native[15, 12] = 0.0 (it is masked, thus zero)

If we used a sub_size of 3, for pixel 0 we we would create a 3x3 sub-array:

        slim[0] = value of first sub-pixel in pixel 0.
        slim[1] = value of first sub-pixel in pixel 1.
        slim[2] = value of first sub-pixel in pixel 2.
0 1 2   slim[3] = value of first sub-pixel in pixel 3.
3 4 5   slim[4] = value of first sub-pixel in pixel 4.
6 7 8   slim[5] = value of first sub-pixel in pixel 5.
        slim[6] = value of first sub-pixel in pixel 6.
        slim[7] = value of first sub-pixel in pixel 7.
        slim[8] = value of first sub-pixel in pixel 8.
Parameters:
  • mask – The ndarray of shape [total_y_pixels, total_x_pixels] containing the bool’s representing the mask, where False signifies an entry is unmasked and used in calculations.

  • pixel_scales – The (y,x) scaled units to pixel units conversion factors of every pixel. If this is input as a float, it is converted to a (float, float) structure.

  • origin – The (y,x) scaled units origin of the mask’s coordinate system.

Methods

all

all_false

Create a mask where all pixels are False and therefore unmasked.

astype

circular

Returns a Mask2D (see Mask2D.__new__) where all False entries are within a circle of input radius.

circular_annular

Returns a Mask2D (see Mask2D.__new__) where all False entries are within an annulus of input inner radius and outer radius.

circular_anti_annular

Returns a Mask2D (see Mask2D.__new__) where all False entries are within an inner circle and second outer circle, forming an inverse annulus.

copy

elliptical

Returns a Mask2D (see Mask2D.__new__) where all False entries are within an ellipse.

elliptical_annular

Returns a Mask2D (see Mask2D.__new__) where all False entries are within an elliptical annulus of input inner and outer scaled major-axis and centre.

flip_hdu_for_ds9

from_cosmic_ray_map_buffed

Returns the mask used for CTI Calibration, which is all False unless specific regions are input for masking.

from_fits

Loads the image from a .fits file.

from_masked_regions

from_pixel_coordinates

Returns a Mask2D (see Mask2D.__new__) where all False entries are defined from an input list of list of pixel coordinates.

from_primary_hdu

Returns an Mask2D by from a PrimaryHDU object which has been loaded via astropy.fits

instance_flatten

Flatten an instance of an autoarray class into a tuple of its attributes (i.e.

instance_unflatten

Unflatten a tuple of attributes (i.e.

invert

manual

Create a Mask2D (see Mask2D.__new__) by inputting the array values in 2D, for example:

mask_new_sub_size_from

Returns the mask on the same scaled coordinate system but with a sub-grid of an inputsub_size.

masked_fpr_and_eper_from

masked_parallel_eper_from

masked_parallel_fpr_from

masked_readout_persistence_from

Read noise persistence is a feature of CCDs whereby the signal from a high signal pixel (e.g.

masked_serial_eper_from

masked_serial_fpr_from

max

min

output_to_fits

Write the 2D Mask to a .fits file.

reshape

sqrt

sum

trimmed_array_from

Map a padded 1D array of values to its original 2D array, trimming all edge values.

unmasked_blurred_array_from

For a padded grid and psf, compute an unmasked blurred image from an unmasked unblurred image.

with_new_array

Copy this object but give it a new array.

Attributes

array

circular_radius

A property that is only computed once per instance and then replaces itself with an ordinary attribute.

derive_grid

derive_indexes

derive_mask

dimensions

dtype

geometry

Return the 2D geometry of the mask, representing its uniform rectangular grid of (y,x) coordinates defined by its shape_native.

hdu_for_output

The mask as a HDU object, which can be output to a .fits file.

imag

is_all_false

Returns False if all pixels in a mask are False, else returns True.

is_all_true

Returns True if all pixels in a mask are True, else returns False.

is_circular

Returns whether the mask is circular or not.

mask

mask_centre

native

Returns the data structure in its native format which contains all unmaksed values to the native dimensions.

ndim

pixel_scale

For a mask with dimensions two or above check that are pixel scales are the same, and if so return this single value as a float.

pixel_scale_header

Returns the pixel scale of the mask as a header dictionary, which can be written to a .fits file.

pixel_scales

pixels_in_mask

The total number of unmasked pixels (values are False) in the mask.

real

shape

shape_native

shape_native_masked_pixels

The (y,x) shape corresponding to the extent of unmasked pixels that go vertically and horizontally across the mask.

shape_slim

The 1D shape of the mask, which is equivalent to the total number of unmasked pixels in the mask.

size

sub_fraction

The fraction of the area of a pixel every sub-pixel contains.

sub_length

The total number of sub-pixels in a give pixel,

sub_pixels_in_mask

The total number of unmasked sub-pixels (values are False) in the mask.

sub_shape_native

sub_shape_slim

The 1D shape of the mask's sub-grid, which is equivalent to the total number of unmasked pixels in the mask.

zoom_centre

zoom_mask_unmasked

The scaled-grid of (y,x) coordinates of every pixel.

zoom_offset_pixels

zoom_offset_scaled

zoom_region

The zoomed rectangular region corresponding to the square encompassing all unmasked values.

zoom_shape_native

classmethod manual(mask, pixel_scales, origin=(0.0, 0.0), invert=False)[source]#

Create a Mask2D (see Mask2D.__new__) by inputting the array values in 2D, for example:

mask=np.array([[False, False],

[True, False]])

mask=[[False, False],

[True, False]]

Parameters:
  • mask – The bool values of the mask input as an ndarray of shape [total_y_pixels, total_x_pixels ]or a list of lists.

  • pixel_scales – The (y,x) arcsecond-to-pixel units conversion factor of every pixel. If this is input as a float, it is converted to a (float, float).

  • origin ((float, float)) – The (y,x) scaled units origin of the mask’s coordinate system.

  • invert – If True, the bool’s of the input mask are inverted, for example False’s become True and visa versa.

classmethod all_false(shape_native, pixel_scales, origin=(0.0, 0.0), invert=False)[source]#

Create a mask where all pixels are False and therefore unmasked.

Parameters:
  • mask – The bool values of the mask input as an ndarray of shape [total_y_pixels, total_x_pixels ]or a list of lists.

  • pixel_scales – The (y,x) arcsecond-to-pixel units conversion factor of every pixel. If this is input as a float, it is converted to a (float, float).

  • origin ((float, float)) – The (y,x) scaled units origin of the mask’s coordinate system.

  • invert – If True, the bool’s of the input mask are inverted, for example False’s become True and visa versa.

classmethod from_cosmic_ray_map_buffed(cosmic_ray_map, settings=<autocti.mask.mask_2d.SettingsMask2D object>)[source]#

Returns the mask used for CTI Calibration, which is all False unless specific regions are input for masking.

Parameters:
  • cosmic_ray_map (array_2d.Array2D) – 2D arrays flagging where cosmic rays on the image.

  • cosmic_ray_parallel_buffer – The number of pixels from each ray pixels are masked in the parallel direction.

  • cosmic_ray_serial_buffer – The number of pixels from each ray pixels are masked in the serial direction.

  • cosmic_ray_diagonal_buffer – The number of pixels from each ray pixels are masked in the digonal up from the parallel + serial direction.

classmethod from_fits(file_path, pixel_scales, hdu=0, origin=(0.0, 0.0), resized_mask_shape=None)[source]#

Loads the image from a .fits file.

Parameters:
  • file_path – The full path of the fits file.

  • hdu – The HDU number in the fits file containing the image image.

  • (float (pixel_scales or) – The arc-second to pixel conversion factor of each pixel.

  • float) – The arc-second to pixel conversion factor of each pixel.

classmethod masked_readout_persistence_from(layout: Layout2D, row_value_list: List[float], readout_persistence_threshold: float, settings: SettingsMask2D, pixel_scales: Union[Tuple[float], Tuple[float, float], float], invert: bool = False) Mask2D[source]#

Read noise persistence is a feature of CCDs whereby the signal from a high signal pixel (e.g. cosmic ray) can persist into the signal of subsequent rows of pixels.

This leads to a ‘streak’ of signal values in the x direction, which typically need to be masked out.

This function produces a read noise persistence mask from a list of row values, where the values are the average signal in each row of the image after other features (e.g. the charge injection) have been removed.

All rows with a signal above an input readout_persistence_threshold are masked out, where this threshold should be estimated from the data itself or based on the CCD’s properties.

Parameters:
  • layout – The layout of the CCD (where the parallel overscan begins and ends, where the charge injection regions are, etc.).

  • row_value_list – The average signal in each row of the image after other features (e.g. the charge injection) have been removed.

  • readout_persistence_threshold – The threshold above which a row is masked out, assuming that this threshold means that a signal is so bright that it must be due to read noise persistence.

  • settings – The settings of the mask (e.g. the number of pixels to mask out).

  • pixel_scales – The pixel scales of the CCD in arc-seconds per pixel, which is passed to the mask.

  • invert – If True, the mask is inverted such that all pixels that are masked are unmasked and visa versa.

Return type:

The read noise persistence mask.