enpt.processors.dead_pixel_correction package

Submodules

enpt.processors.dead_pixel_correction.dead_pixel_correction module

EnPT ‘dead pixel correction’ module.

Performs the Dead Pixel Correction using a linear interpolation in spectral dimension.

class enpt.processors.dead_pixel_correction.dead_pixel_correction.Dead_Pixel_Corrector(algorithm: str = 'spatial', interp_spectral: str = 'linear', interp_spatial: str = 'linear', CPUs: int = None, logger: Logger = None)[source]

Bases: object

EnPT Dead Pixel Correction class.

The EnPT dead pixel correction uses the pixel masks provided by DLR and interpolates the EnMAP image data at the indicated dead pixel positions. It supports two interpolation algorithms:

  1. spectral interpolation
    • Interpolates the data in the spectral domain.

    • Points outside the data range are extrapolated.

  2. spatial interpolation
    • Interpolates the data spatially.

    • Remaining missing data positions (e.g., outermost columns) are spectrally interpolated.

Get an instance of Dead_Pixel_Corrector.

Parameters:
  • algorithm – algorithm how to correct dead pixels ‘spectral’: interpolate in the spectral domain ‘spatial’: interpolate in the spatial domain

  • interp_spectral – spectral interpolation algorithm (‘linear’, ‘quadratic’, ‘cubic’)

  • interp_spatial – spatial interpolation algorithm (‘linear’, ‘bilinear’, ‘cubic’, ‘spline’)

  • CPUs – number of CPUs to use for interpolation (only relevant if algorithm = ‘spatial’)

  • logger

_interpolate_nodata_spatially(image2correct: GeoArray, deadpixel_map: GeoArray)[source]
_interpolate_nodata_spectrally(image2correct: GeoArray, deadpixel_map: GeoArray)[source]
static _validate_inputs(image2correct: GeoArray, deadpixel_map: GeoArray)[source]
correct(image2correct: ndarray | GeoArray, deadpixel_map: ndarray | GeoArray)[source]

Run the dead pixel correction.

Parameters:
  • image2correct – image to correct

  • deadpixel_map – dead pixel map (2D (bands x columns) or 3D (rows x columns x bands)

Returns:

corrected image

enpt.processors.dead_pixel_correction.dead_pixel_correction._get_baddata_mask(data: ndarray, nodata: ndarray | Number = nan, transpose_inNodata: bool = False)[source]
enpt.processors.dead_pixel_correction.dead_pixel_correction.interp_nodata_along_axis(data, axis=0, nodata: ndarray | Number = nan, method: str = 'linear', **kw)[source]

Interpolate a 2D or 3D array along the given axis (based on scipy.interpolate.make_interp_spline).

Parameters:
  • data – data to interpolate

  • axis – axis to interpolate (0: along columns; 1: along rows, 2: along bands)

  • nodata – nodata array in the shape of data or nodata value

  • method – interpolation method (‘linear’, ‘quadratic’, ‘cubic’)

  • kw – keyword arguments to be passed to scipy.interpolate.make_interp_spline

Returns:

interpolated array

enpt.processors.dead_pixel_correction.dead_pixel_correction.interp_nodata_along_axis_2d(data_2d: ndarray, axis: int = 0, nodata: ndarray | Number = nan, method: str = 'linear', **kw)[source]

Interpolate a 2D array along the given axis (based on scipy.interpolate.make_interp_spline).

Parameters:
  • data_2d – data to interpolate

  • axis – axis to interpolate (0: along columns; 1: along rows)

  • nodata – nodata array in the shape of data or nodata value

  • method – interpolation method (‘linear’, ‘quadratic’, ‘cubic’)

  • kw – keyword arguments to be passed to scipy.interpolate.make_interp_spline

Returns:

interpolated array

enpt.processors.dead_pixel_correction.dead_pixel_correction.interp_nodata_spatially_2d(data_2d: ndarray, axis: int = 0, nodata: ndarray | Number = nan, method: str = 'linear', fill_value: float = nan, implementation: str = 'pandas') ndarray[source]

Interpolate a 2D array spatially.

NOTE: If data_2d contains NaN values that are unlabelled by a given nodata array,

they are also overwritten in the pandas implementation.

Parameters:
  • data_2d – data to interpolate

  • axis – axis to interpolate (0: along columns; 1: along rows)

  • nodata – nodata array in the shape of data or nodata value

  • method – interpolation method - if implementation=’scipy’: ‘linear’, ‘nearest’, ‘cubic’ - if implementation=’pandas’: ‘linear’, ‘nearest’, ‘slinear’, ‘quadratic’, ‘cubic’, etc.

  • fill_value – value to fill into positions where no interpolation is possible

  • implementation – ‘scipy’: interpolation based on scipy.interpolate.griddata ‘pandas’: interpolation based on pandas.core.resample.Resampler.interpolate

Returns:

interpolated array

enpt.processors.dead_pixel_correction.dead_pixel_correction.interp_nodata_spatially_3d(data_3d: ndarray, axis: int = 0, nodata: ndarray | Number = nan, method: str = 'linear', fill_value: float = nan, implementation: str = 'pandas', CPUs: int = None) ndarray[source]

Interpolate a 3D array spatially, band-for-band.

Parameters:
  • data_3d – data to interpolate

  • axis – axis to interpolate (0: along columns; 1: along rows)

  • nodata – nodata array in the shape of data or nodata value

  • method – interpolation method - if implementation=’scipy’: ‘linear’, ‘nearest’, ‘cubic’ - if implementation=’pandas’: ‘linear’, ‘nearest’, ‘slinear’, ‘quadratic’, ‘cubic’, etc.

  • fill_value – value to fill into positions where no interpolation is possible

  • implementation – ‘scipy’: interpolation based on scipy.interpolate.griddata ‘pandas’: interpolation based on pandas.core.resample.Resampler.interpolate

  • CPUs – number of CPUs to use

Returns:

interpolated array

Module contents

EnPT ‘dead pixel correction’ module for detecting and correcting dead pixels.