vip_hci.preproc package

Submodules

vip_hci.preproc.badframes module

Module with functions for outlier frame detection.

vip_hci.preproc.badframes.cube_detect_badfr_correlation(array, frame_ref, crop_size=30, dist='pearson', percentile=20, threshold=None, mode='full', inradius=None, width=None, plot=True, verbose=True, full_output=False)[source]

Returns the list of bad frames from a cube by measuring the distance (similarity) or correlation of the frames (cropped to a 30x30 subframe) wrt a reference frame from the same cube. Then the distance/correlation level is thresholded (percentile parameter) to find the outliers. Should be applied on a recentered cube.

Parameters:
  • array (numpy ndarray) – Input 3d array, cube.

  • frame_ref (int or 2d array) – Index of the frame that will be used as a reference or 2d reference array.

  • crop_size (int, optional) – Size in pixels of the square subframe to be analyzed.

  • dist ({'sad','euclidean','mse','pearson','spearman','ssim'}, str optional) – One of the similarity or dissimilarity measures from function vip_hci.stats.distances.cube_distance.

  • percentile (int, optional) – The percentage of frames that will be discarded, if threshold is not provided.

  • threshold (None or float, optional) – If provided, corresponds to the threshold ‘distance’ value above/below which (depending on index of similarity/dissimilarity resp.) will be discarded. If not None, supersedes ‘percentile’.

  • mode ({'full','annulus'}, string optional) – Whether to use the full frames or a centered annulus.

  • inradius (None or int, optional) – The inner radius when mode is ‘annulus’.

  • width (None or int, optional) – The width when mode is ‘annulus’.

  • plot (bool, optional) – If true it plots the mean fluctuation as a function of the frames and the boundaries.

  • verbose (bool, optional) – Whether to print to stdout or not.

  • full_output (bool, optional) – Whether to also return the array of distances.

Returns:

  • good_index_list (numpy ndarray) – 1d array of good indices.

  • bad_index_list (numpy ndarray) – 1d array of bad frames indices.

  • distances (numpy ndarray) – [if full_output=True] 1d array with the measured distances

vip_hci.preproc.badframes.cube_detect_badfr_ellipticity(array, fwhm, crop_size=30, roundlo=-0.2, roundhi=0.2, plot=True, verbose=True)[source]

Returns the list of bad frames from a cube by measuring the PSF ellipticity of the central source. Should be applied on a recentered cube.

Parameters:
  • array (numpy ndarray) – Input 3d array, cube.

  • fwhm (float) – FWHM size in pixels.

  • crop_size (int, optional) – Size in pixels of the square subframe to be analyzed.

  • roundlo (float, optional) – Lower and higher bounds for the ellipticity. See Notes below for details.

  • roundhi (float, optional) – Lower and higher bounds for the ellipticity. See Notes below for details.

  • plot (bool, optional) – If true it plots the central PSF roundness for each frame.

  • verbose (bool, optional) – Whether to print to stdout or not.

Returns:

  • good_index_list (numpy ndarray) – 1d array of good indices.

  • bad_index_list (numpy ndarray) – 1d array of bad frames indices.

Note

From photutils.DAOStarFinder documentation: DAOFIND calculates the object roundness using two methods. The ‘roundlo’ and ‘roundhi’ bounds are applied to both measures of roundness. The first method (‘roundness1’; called ‘SROUND’ in DAOFIND) is based on the source symmetry and is the ratio of a measure of the object’s bilateral (2-fold) to four-fold symmetry. The second roundness statistic (‘roundness2’; called ‘GROUND’ in DAOFIND) measures the ratio of the difference in the height of the best fitting Gaussian function in x minus the best fitting Gaussian function in y, divided by the average of the best fitting Gaussian functions in x and y. A circular source will have a zero roundness. A source extended in x or y will have a negative or positive roundness, respectively.

vip_hci.preproc.badframes.cube_detect_badfr_pxstats(array, mode='annulus', in_radius=10, width=10, top_sigma=1.0, low_sigma=1.0, window=None, plot=True, verbose=True)[source]

Returns the list of bad frames from a cube using the px statistics in a centered annulus or circular aperture. Frames that are more than a few standard deviations discrepant are rejected. Should be applied on a recentered cube.

Parameters:
  • array (numpy ndarray) – Input 3d array, cube.

  • mode ({'annulus', 'circle'}, string optional) – Whether to take the statistics from a circle or an annulus.

  • in_radius (int optional) – If mode is ‘annulus’ then ‘in_radius’ is the inner radius of the annular region. If mode is ‘circle’ then ‘in_radius’ is the radius of the aperture.

  • width (int optional) – Size of the annulus. Ignored if mode is ‘circle’.

  • top_sigma (int, optional) – Top boundary for rejection.

  • low_sigma (int, optional) – Lower boundary for rejection.

  • window (int, optional) – Window for smoothing the mean and getting the rejection statistic. If None, it is defined as n_frames//3.

  • plot (bool, optional) – If true it plots the mean fluctuation as a function of the frames and the boundaries.

  • verbose (bool, optional) – Whether to print to stdout or not.

Returns:

  • good_index_list (numpy ndarray) – 1d array of good indices.

  • bad_index_list (numpy ndarray) – 1d array of bad frames indices.

vip_hci.preproc.badpixremoval module

Module with functions for correcting bad pixels in cubes.

[AAC01] (1,2,3,4,5)
Aach & Metzler 2001
Defect interpolation in digital radiography how object-oriented transform coding helps
SPIE, Proceedings Volume 4322, Medical Imaging 2001: Image Processing
vip_hci.preproc.badpixremoval.cube_fix_badpix_annuli(array, fwhm, cy=None, cx=None, sig=5.0, bpm_mask=None, protect_mask=0, excl_mask=None, r_in_std=50, r_out_std=None, verbose=True, half_res_y=False, min_thr=None, max_thr=None, min_thr_np=None, bad_values=None, full_output=False)[source]

Function to correct the bad pixels annulus per annulus (centered on the provided location of the star), in an input frame or cube. This function is faster than cube_fix_badpix_clump; hence to be preferred in all cases where there is only one bright source with circularly symmetric PSF. The bad pixel values are replaced by: ann_median + random_poisson; where ann_median is the median of the annulus, and random_poisson is random noise picked from a Poisson distribution centered on ann_median.

Parameters:
  • array (3D or 2D array) – Input 3d cube or 2d image.

  • fwhm (float or 1D array) – Vector containing the full width half maximum of the PSF in pixels, for each channel (cube_like); or single value (frame_like)

  • cy (None, float or 1D array, optional) – If None: will use the barycentre of the image found by photutils.centroid_com() If floats: coordinates of the center, assumed to be the same in all frames if the input is a cube. If 1D arrays: they must be the same length as the 0th dimension of the input cube.

  • cx (None, float or 1D array, optional) – If None: will use the barycentre of the image found by photutils.centroid_com() If floats: coordinates of the center, assumed to be the same in all frames if the input is a cube. If 1D arrays: they must be the same length as the 0th dimension of the input cube.

  • sig (Float scalar, optional) – Number of stddev above or below the median of the pixels in the same annulus, to consider a pixel as bad.

  • bpm_mask (3D or 2D array, opt) – Input bad pixel array. If 2D and array is 3D: should have same last 2 dimensions as array. If 3D, should have exact same dimensions as input array. If not provided, the algorithm will attempt to identify bad pixel clumps automatically.

  • protect_mask (int or float, optional) – If larger than 0, radius of a circular aperture (at the center of the frames) in which no bad pixels will be identified. This can be useful to protect the star and vicinity.

  • excl_mask (numpy ndarray, optional) – Binary mask with 1 in areas that should not be considered as good neighbouring pixels during the identification of bad pixels. These should not be considered as bad pixels to be corrected neither (i.e. different to bpm_mask).

  • r_in_std (float or None, optional) – Inner radius (in pixels) of the annulus used for the calculation of the standard deviation of the background noise - used as a min threshold when identifying bad pixels. Default: 50.

  • r_out_std (float or None, optional) – Outer radius in pixels of the annulus used for the calculation of the standard deviation of the background noise - used as a min threshold when identifying bad pixels. If set to None, the default will be to consider the largest annulus starting at r_in_std which fits within the frame.

  • verbose (bool, {False, True}, optional) – Whether to print out the number of bad pixels in each frame.

  • half_res_y (bool, {True,False}, optional) – Whether the input data have only half the angular resolution vertically compared to horizontally (e.g. SINFONI data). The algorithm will then correct the bad pixels every other row.

  • min_thr ({None,float}, optional) – Any pixel whose value is lower (resp. larger) than this threshold will be automatically considered bad and hence sigma_filtered. If None, it is not used.

  • max_thr ({None,float}, optional) – Any pixel whose value is lower (resp. larger) than this threshold will be automatically considered bad and hence sigma_filtered. If None, it is not used.

  • min_thr_np ({None, float}, optional) – Any pixel whose value is lower than this threshold will be automatically considered bad and hence sigma_filtered, EVEN if located within the radius of protect_mask.

  • bad_values (list or None, optional) – If not None, should correspond to a list of known bad values (e.g. 0). These pixels will be added to the input bad pixel map.

  • full_output (bool, {False,True}, optional) – Whether to return as well the cube of bad pixel maps and the cube of defined annuli.

Returns:

  • array_corr (2d or 3d array) – The bad pixel corrected frame/cube.

  • bpix_map (2d or 3d array) – [full_output=True] The bad pixel map or the cube of bpix maps

  • ann_frame_cumul (2 or 3d array) – [full_output=True] The cube of defined annuli

vip_hci.preproc.badpixremoval.cube_fix_badpix_clump(array, bpm_mask=None, correct_only=False, cy=None, cx=None, fwhm=4.0, sig=4.0, protect_mask=0, excl_mask=None, half_res_y=False, min_thr=None, max_nit=15, mad=True, bad_values=None, verbose=True, full_output=False, nproc=1)[source]

Function to identify and correct clumps of bad pixels. Very fast when a bad pixel map is provided. If a bad pixel map is not provided, the bad pixel clumps will be searched iteratively and replaced by the median of good neighbouring pixel values, when enough of them are available. The size of the box is set by the closest odd integer larger than fwhm (to avoid accidentally replacing point sources).

Parameters:
  • array (3D or 2D array) – Input 3d cube or 2d image.

  • bpm_mask (3D or 2D array, opt) – Input bad pixel array. If 2D and array is 3D: should have same last 2 dimensions as array. If 3D, should have exact same dimensions as input array. If not provided, the algorithm will attempt to identify bad pixel clumps automatically.

  • correct_only (bool, opt) – If True and bpix_map is provided, will only correct for provided bad pixels. Else, the algorithm will determine (more) bad pixels.

  • cy (float or 1D array. opt) – Vector with approximate y and x coordinates of the star for each channel (cube_like), or single 2-elements vector (frame_like). These will be used if bpix_map is None and protect_mask set to True. If left to None, default values will correspond to the central pixel coordinates.

  • cx (float or 1D array. opt) – Vector with approximate y and x coordinates of the star for each channel (cube_like), or single 2-elements vector (frame_like). These will be used if bpix_map is None and protect_mask set to True. If left to None, default values will correspond to the central pixel coordinates.

  • fwhm (float or 1D array, opt) – Vector containing the full width half maximum of the PSF in pixels, for each channel (cube_like); or single value (frame_like). Should be provided if bpix map is None.

  • sig (float, optional) – Value representing the number of “sigmas” above or below the “median” of the neighbouring pixel, to consider a pixel as bad. See details on parameter “m” of function reject_outlier.

  • protect_mask (int or float, optional) – If larger than 0, radius of a circular aperture (at the center of the frames) in which no bad pixels will be identified. This can be useful to protect the star and vicinity.

  • excl_mask (numpy ndarray, optional) – Binary mask with same dimensions as array, with 1 in areas that should not be considered as good neighbouring pixels during the identification of bad pixels. These should not be considered as bad pixels to be corrected neither (i.e. different to bpm_mask).

  • half_res_y (bool, {True,False}, optional) – Whether the input data has only half the angular resolution vertically compared to horizontally (e.g. the case of SINFONI data); in other words there are always 2 rows of pixels with exactly the same values. The algorithm will just consider every other row (hence making it twice faster), then apply the bad pixel correction on all rows.

  • min_thr (float, tuple or None, opt) – If a float is provided, corresponds to a minimum absolute threshold below which pixels are not considered bad (can be used to avoid the identification of bad pixels within noise). If a tuple of 2 values, corresponds to the range of values within which not to consider a pixel as bad. (e.g. (-0.1, 10.)).

  • max_nit (float, optional) – Maximum number of iterations on a frame to correct bpix. Typically, it should be set to less than ny/2 or nx/2. This is a mean of precaution in case the algorithm gets stuck with 2 neighbouring pixels considered bpix alternately on two consecutively iterations hence leading to an infinite loop (very very rare case).

  • mad ({False, True}, bool optional) – If True, the median absolute deviation will be used instead of the standard deviation.

  • bad_values (list or None, optional) – If not None, should correspond to a list of known bad values (e.g. 0). These pixels will be added to the input bad pixel map.

  • verbose (bool, {False,True}, optional) – Whether to print the number of bad pixels and number of iterations required for each frame.

  • full_output (bool, {False,True}, optional) – Whether to return as well the cube of bad pixel maps and the cube of defined annuli.

  • nproc (int, optional) – This feature is added following ADACS update. Refers to the number of processors available for calculations. Choosing a number >1 enables multiprocessing for the correction of frames.

Returns:

  • array_corr (2d or 3d array) – The bad pixel corrected frame/cube.

  • bpix_map (2d or 3d array) – [full_output=True] The bad pixel map or the cube of bpix maps

vip_hci.preproc.badpixremoval.cube_fix_badpix_ifs(array, lbdas, fluxes=None, mask=None, cy=None, cx=None, clumps=True, sigma_clip=3, num_neig=5, size=5, protect_mask=0, mad=False, fwhm=4, min_thr=None, max_nit=15, imlib='vip-fft', interpolation='lanczos4', ignore_nan=True, verbose=True, full_output=False)[source]

Function to identify and correct bad pixels in an IFS cube, leveraging on the radial expansion of the PSF with wavelength. Bad pixel identification is done with either the cube_fix_badpix_isolated or the cube_fix_badpix_clump function in PSF subtracted frames (through SDI).

Parameters:
  • array (3D or 4D array) – Input 3D (spectral) or 4D (spectral+temporal) cube. In the latter case, dimensions should be spectral x temporal x vertical x horizontal.

  • lbdas (1d array or list) – Vector with the wavelengths, used for first guess on scaling factor.

  • fluxes (1d array or list, optional) – Vector with the (unsaturated) fluxes at the different wavelengths, used for first guess on flux factor.

  • mask (2D-array, opt) – Binary mask, with ones where the residual intensities should be evaluated. If None is provided, the whole field is used.

  • cy (None, float or 1D array, optional) – If None: will use the barycentre of the image found by photutils.centroid_com() If floats: coordinates of the center, assumed to be the same in all frames if the input is a cube. If 1D arrays: they must be the same length as the 0th dimension of the input cube.

  • cx (None, float or 1D array, optional) – If None: will use the barycentre of the image found by photutils.centroid_com() If floats: coordinates of the center, assumed to be the same in all frames if the input is a cube. If 1D arrays: they must be the same length as the 0th dimension of the input cube.

  • clumps (bool, optional) – Whether to use cube_fix_badpix_clump (True) or cube_fix_badpix_isolated (False) in the SDI residual cube.

  • sigma_clip (int, optional) – In case no bad pixel mask is provided all the pixels above and below sigma_clip*STDDEV will be marked as bad.

  • num_neig (int, optional) – The side of the square window around each pixel where the sigma clipped statistics are calculated (STDDEV and MEDIAN). If the value is equal to 0 then the statistics are computed in the whole frame.

  • size (odd int, optional) – The size the box (size x size) of adjacent pixels for the median filter.

  • protect_mask (int or float, optional) – If larger than 0, radius of a circular aperture (at the center of the frames) in which no bad pixels will be identified. This can be useful to protect the star and vicinity.

  • mad ({False, True}, bool optional) – If True, the median absolute deviation will be used instead of the standard deviation.

  • fwhm (float or 1D array, opt) – Vector containing the full width half maximum of the PSF in pixels, for each channel (cube_like); or single value (frame_like). Shouod be provided if bpix map is None.

  • min_thr (float, tuple or None, opt) – If a float is provided, corresponds to a minimum absolute threshold below which pixels are not considered bad in the residua images (can be used to avoid the identification of bad pixels within noise). If a tuple of 2 values, corresponds to the range of values within which not to consider a pixel as bad, in the residual images (e.g. (-1, 5)).

  • max_nit (float, optional) – Maximum number of iterations on a frame to correct bpix. Typically, it should be set to less than ny/2 or nx/2. This is a mean of precaution in case the algorithm gets stuck with 2 neighbouring pixels considered bpix alternately on two consecutively iterations hence leading to an infinite loop (very very rare case).

  • ignore_nan (bool, optional) – Whether to not consider NaN values as bad pixels. If False, will also correct them.

  • verbose (bool, optional) – If True additional information will be printed.

  • full_output (bool, {False,True}, optional) – Whether to return as well the cube of bad pixel maps and the cube of defined annuli.

Returns:

  • array_out (numpy ndarray) – Cube with bad pixels corrected.

  • bpm_mask (2d or 3d array [if full_output is True]) – The bad pixel map or the cube of bpix maps

  • array_res (2d or 3d array [if full_output is True]) – SDI-residual cube in which bad pixels are identified

vip_hci.preproc.badpixremoval.cube_fix_badpix_interp(array, bpm_mask, mode='fft', excl_mask=None, fwhm=4.0, kernel_sz=None, psf=None, half_res_y=False, nit=500, tol=1, nproc=1, full_output=False, **kwargs)[source]

Function to correct clumps of bad pixels by interpolation with either a user-defined kernel (through astropy.convolution) or through the FFT-based algorithm described in [AAC01]. A bad pixel map must be provided (e.g. found with function cube_fix_badpix_clump).

Parameters:
  • array (3D or 2D array) – Input 3d cube or 2d image.

  • bpm_mask (3D or 2D array) – Input bad pixel array. Should have same x,y dimensions as array. If 2D, but input array is 3D, the same bpix_map will be assumed for all frames.

  • mode (str, optional {'fft', 'gauss', 'psf'}) – Can be either a 2D Gaussian (‘gauss’) or an input normalized PSF (‘psf’).

  • excl_mask (3D or 2D array, optional) – [Only used if mode != ‘fft’] Input exclusion mask array. Pixels in the exclusion mask will neither be used for interpolation, nor replaced as bad pixels. excl_mask should have same x,y dimensions as array. If 2D, but input array is 3D, the same exclusion mask will be assumed for all frames.

  • fwhm (float, 1D array or tuple of 2 floats, opt) – If mode is ‘gauss’, the fwhm of the Gaussian.

  • kernel_sz (int or None, optional) – Size of the kernel in pixels for 2D Gaussian and Moffat convolutions. If None, astropy.convolution will automatically consider 8*fwhm kernel sizes.

  • psf (2D or 3D array, optional) – If mode is ‘psf’, a normalized PSF array. If a 3D cube is provided (e.g. for spectral cubes), the first dimension should match that of the input array (which should also be 3D). Else, the same 2D PSF kernel will be for all input frames, whether the input is 2D or 3D. If half_res_y is True, psf should be provided vertically squashed.

  • half_res_y (bool, {True,False}, optional) – Whether the input data has only half the angular resolution vertically compared to horizontally (e.g. the case for some IFUs); in other words there are always 2 rows of pixels with exactly the same values. If so, the Gaussian kernel will also be squashed vertically by a factor 2.

  • nit (int or list of int, optional) – For FFT-based iterative interpolation, the number of iterations to use. If a list is provided, a list of bad pixel corrected frames/cubes is returned.

  • tol (float) – Tolerance in terms of E_g (see [AAC01]). The iterative process is stopped if the error E_g gets lower than this tolerance.

  • nproc (None or int, optional) – Number of processes for parallel computing. If None the number of processes will be set to (cpu_count()/2). Note: only used for input 3D cube and ‘fft’ mode.

  • full_output (bool, {False,True}, optional) – In the case of FT-based interpolation, whether to return as well the reconstructed images.

  • **kwargs (dict) – Passed through to the astropy.convolution.convolve or convolve_fft function.

Returns:

  • array_corr (2d or 3d array;) – The bad pixel corrected frame/cube.

  • recon_cube (2d or 3d array;) – [full_output=True & mode=’fft’] The reconstructed frame/cube. If nit is a list, a list of reconstructed frames/cubes is returned.

vip_hci.preproc.badpixremoval.cube_fix_badpix_isolated(array, bpm_mask=None, correct_only=False, sigma_clip=3, num_neig=5, size=5, frame_by_frame=False, protect_mask=0, excl_mask=None, cxy=None, mad=False, ignore_nan=True, verbose=True, full_output=False, nproc=1)[source]

Corrects the bad pixels, marked in the bad pixel mask. The bad pixel is replaced by the median of the adjacent pixels. This function is very fast but works only with isolated (sparse) pixels.

Parameters:
  • array (numpy ndarray) – Input 3d array.

  • bpm_mask (numpy ndarray, optional) – Input bad pixel map. Zeros frame where the bad pixels have a value of 1. If None is provided a bad pixel map will be created per frame using sigma clip statistics.

  • correct_only (bool, opt) – If True and bpix_map is provided, will only correct for provided bad pixels. Else, the algorithm will determine (more) bad pixels.

  • sigma_clip (int, optional) – In case no bad pixel mask is provided all the pixels above and below sigma_clip*STDDEV will be marked as bad.

  • num_neig (int, optional) – The side of the square window around each pixel where the sigma clipped statistics are calculated (STDDEV and MEDIAN). If the value is equal to 0 then the statistics are computed in the whole frame.

  • size (odd int, optional) – The size the box (size x size) of adjacent pixels for the median filter.

  • frame_by_frame (bool, optional) – Whether to correct bad pixels frame by frame in the cube. By default it is set to False; the bad pixels are computed on the mean frame of the stack (faster but not necessarily optimal).

  • protect_mask (int or float, optional) – If larger than 0, radius of a circular aperture (at the center of the frames) in which no bad pixels will be identified. This can be useful to protect the star and vicinity.

  • excl_mask (numpy ndarray, optional) – Binary mask with 1 in areas that should not be considered as good neighbouring pixels during the identification of bad pixels. These should not be considered as bad pixels to be corrected neither (i.e. different to bpm_mask).

  • cxy (None, tuple or 2d numpy ndarray) – If protect_mask is True, this is the location of the star centroid in the images. If None, assumes the star is already centered. If a tuple, the location of the star is assumed to be the same in all frames of the cube. If a (n_frames x 2) ndarray, it should contain the xy location of the star in each frame.

  • mad ({False, True}, bool optional) – If True, the median absolute deviation will be used instead of the standard deviation.

  • ignore_nan (bool, optional) – Whether to not consider NaN values as bad pixels. If False, will also correct them.

  • verbose (bool, optional) – If True additional information will be printed.

  • full_output (bool, {False,True}, optional) – Whether to return as well the cube of bad pixel maps and the cube of defined annuli.

  • nproc (int, optional) – This feature is added following ADACS update. Refers to the number of processors available for calculations. Choosing a number >1 enables multiprocessing for the correction of frames. This happens only when ``frame_by_frame=True’’.

Returns:

  • array_out (numpy ndarray) – Cube with bad pixels corrected.

  • bpm_mask (2d or 3d array [if full_output is True]) – The bad pixel map or the cube of bad pixel maps

vip_hci.preproc.badpixremoval.frame_fix_badpix_fft(array, bpm_mask, nit=500, tol=1, pad_fac=2, verbose=True, full_output=False)[source]

Function to interpolate bad pixels with the FFT-based algorithm in [AAC01].

Parameters:
  • array (2D ndarray) – Input image.

  • bpm_mask (2D ndarray) – Bad pixel map.

  • nit (int or list of int, optional) – The number of iterations to use. If a list is provided, a list of bad pixel corrected frames/cubes is returned.

  • tol (float, opt) – Tolerance in terms of E_g (see [AAC01]). The iterative process is stopped if the error E_g gets lower than this tolerance.

  • pad_fac (int or float, opt) – Padding factor before calculating 2D-FFT.

  • verbose (bool) – Whether to print additional information during processing, incl. progress bar.

  • full_output (bool) – Whether to also return the reconstructed estimate f_hat of the input array.

Returns:

  • array_corr (2D ndarray or list of 2D ndarray) – Image in which the bad pixels have been interpolated.

  • f_est (2D ndarray or list of 2D ndarray) – [full_output=True] Reconstructed estimate (f_hat in [AAC01]) of the input array

vip_hci.preproc.badpixremoval.frame_fix_badpix_isolated(array, bpm_mask=None, correct_only=False, sigma_clip=3, num_neig=5, size=5, protect_mask=0, excl_mask=None, cxy=None, mad=False, ignore_nan=True, verbose=True, full_output=False)[source]

Corrects the bad pixels, marked in the bad pixel mask. The bad pixel is replaced by the median of the adjacent pixels. This function is very fast but works only with isolated (sparse) pixels.

Parameters:
  • array (numpy ndarray) – Input 2d array.

  • bpm_mask (numpy ndarray, optional) – Input bad pixel map. Should have same size as array. Binary map where the bad pixels have a value of 1. If None is provided a bad pixel map will be created using sigma clip statistics.

  • correct_only (bool, opt) – If True and bpix_map is provided, will only correct for provided bad pixels. Else, the algorithm will determine (more) bad pixels.

  • sigma_clip (int, optional) – In case no bad pixel mask is provided all the pixels above and below sigma_clip*STDDEV will be marked as bad.

  • num_neig (int, optional) – The side of the square window around each pixel where the sigma clipped statistics are calculated (STDDEV and MEDIAN). If the value is equal to 0 then the statistics are computed in the whole frame.

  • size (odd int, optional) – The size the box (size x size) of adjacent pixels for the median filter.

  • protect_mask (int or float, optional) – If larger than 0, radius of a circular aperture (at the center of the frames) in which no bad pixels will be identified. This can be useful to protect the star and vicinity.

  • excl_mask (numpy ndarray, optional) – Binary mask with 1 in areas that should not be considered as good neighbouring pixels during the identification of bad pixels. These should not be considered as bad pixels to be corrected neither (i.e. different to bpm_mask).

  • cxy (None or tuple) – If protect_mask is True, this is the location of the star centroid in the images. If None, assumes the star is already centered. If a tuple, the location of the star is assumed to be the same in all frames of the cube.

  • mad ({False, True}, bool optional) – If True, the median absolute deviation will be used instead of the standard deviation.

  • ignore_nan (bool, optional) – Whether to not consider NaN values as bad pixels. If False, will also correct them.

  • verbose (bool, optional) – If True additional information will be printed.

  • full_output (bool, {False,True}, optional) – Whether to return as well the cube of bad pixel maps and the cube of defined annuli.

Returns:

  • frame (numpy ndarray) – Frame with bad pixels corrected.

  • bpm_mask (2d array) – The bad pixel map

vip_hci.preproc.cosmetics module

Module with cosmetics procedures. Contains the function for bad pixel fixing. Also functions for cropping cubes.

vip_hci.preproc.cosmetics.approx_stellar_position(cube, fwhm, return_test=False, verbose=False)[source]

Finds the approximate coordinates of the star, assuming it is the brightest signal in the images. The algorithm can handle images dominated by noise, since outliers are corrected based on the position of ths star in other channels.

Parameters:
  • obj_tmp (numpy ndarray) – Input 3d cube

  • fwhm (float or array 1D) – Input full width half maximum value of the PSF for each channel. This will be used as the standard deviation for Gaussian kernel of the Gaussian filtering. If float, it is assumed the same for all channels.

  • return_test (bool, optional) – Whether the test result vector (a bool vector with whether the star centroid could be find in the corresponding channel) should be returned as well, along with the approx stellar coordinates.

  • verbose (bool, optional) – Chooses whether to print some additional information.

Returns:

  • star_approx_idx (2d numpy array) – Array of y and x approximate coordinates of the star in each channel of the cube. Dimensions are nz x 2.

  • test_result (1d numpy array) – [return_test=True] It also returns the test result vector.

vip_hci.preproc.cosmetics.cube_correct_nan(cube, neighbor_box=3, min_neighbors=3, verbose=False, half_res_y=False, nproc=1)[source]

Sigma filtering of nan pixels in a whole frame or cube. Tested on SINFONI data.

Parameters:
  • cube (cube_like) – Input 3d or 2d array.

  • neighbor_box (int, optional) – The side of the square window around each pixel where the sigma and median are calculated for the nan pixel correction.

  • min_neighbors (int, optional) – Minimum number of good neighboring pixels to be able to correct the bad/nan pixels.

  • verbose (bool, optional) – Whether to print more information or not during processing

  • half_res_y (bool, optional) – Whether the input data have every couple of 2 rows identical, i.e. there is twice less angular resolution vertically than horizontally (e.g. SINFONI data). The algorithm goes twice faster if this option is rightfully set to True.

  • nproc (None or int, optional) – Number of CPUs for multiprocessing (only used for 3D input). If None, will automatically set it to half the available number of CPUs.

Returns:

obj_tmp – Output cube with corrected nan pixels in each frame

Return type:

numpy ndarray

vip_hci.preproc.cosmetics.cube_crop_frames(array, size, xy=None, force=False, verbose=True, full_output=False)[source]

Crops frames in a cube (3d or 4d array).

Parameters:
  • array (numpy ndarray) – Input 3d or 4d array.

  • size (int) – Size of the desired central sub-array in each frame, in pixels.

  • xy (tuple of ints) – X, Y coordinates of new frame center. If you are getting the coordinates from ds9 subtract 1, python has 0-based indexing.

  • force (bool, optional) – size and the original size of the frames must be both even or odd. With force set to True this condition can be avoided.

  • verbose (bool optional) – If True message of completion is showed.

  • full_output (bool optional) – If true, returns cenx and ceny in addition to array_view.

Returns:

array_out – Cube with cropped frames.

Return type:

numpy ndarray

vip_hci.preproc.cosmetics.cube_drop_frames(array, n, m, parallactic=None, verbose=True)[source]

Slice the cube so that all frames between n and m are kept.

Operates on axis 0 for 3D cubes, and on axis 1 for 4D cubes. This returns a modified copy of array. The indices n and m are included and 1-based.

Parameters:
  • array (numpy ndarray) – Input cube, 3d or 4d.

  • n (int) – 1-based index of the first frame to be kept. Frames before this one are dropped.

  • m (int) – 1-based index of the last frame to be kept. Frames after this one are dropped.

  • parallactic (1d ndarray, optional) – parallactic angles vector. If provided, a modified copy of parallactic is returned additionally.

Returns:

  • array_view (numpy ndarray) – Cube with new size.

  • parallactic (1d numpy ndarray) – [parallactic != None] New parallactic angles.

vip_hci.preproc.cosmetics.frame_crop(array, size, cenxy=None, force=False, verbose=True)[source]

Crops a square frame (2d array). Uses the get_square function.

Parameters:
  • array (numpy ndarray) – Input frame.

  • size (int, odd) – Size of the subframe.

  • cenxy (tuple, optional) – Coordinates of the center of the subframe.

  • force (bool, optional) – Size and the size of the 2d array must be both even or odd. With force set to False, the requested size is flexible (i.e. +1 can be applied to requested crop size for its parity to match the input size). If force set to True, the requested crop size is enforced, even if parities do not match (warnings are raised!).

  • verbose (bool optional) – If True, a message of completion is shown.

Returns:

array_view – Sub array.

Return type:

numpy ndarray

vip_hci.preproc.cosmetics.frame_pad(array, fac, fillwith=0, loc=0, scale=1, keep_parity=True, full_output=False)[source]

Pads a frame (2d array) equally on each sides, where the final frame size is set by a multiplicative factor applied to the original size. The padding is set by fillwith, which can be either a fixed value or white noise, characterized by (loc, scale). Uses the get_square function.

Parameters:
  • array (2D numpy ndarray) – Input frame.

  • fac (float > 1 or tuple of 2 floats > 1.) – Ratio of the size between padded and input frame. If a tuple, corresponds to padding factors along the y and x dimensions resp.

  • fillwith (float or str, optional) – If a float or np.nan: value used for padding. If str, must be ‘noise’, which will inject white noise, using loc and scale parameters.

  • loc (float, optional) – If padding noise, mean of the white noise.

  • scale (float, optional) – If padding noise, standard deviation of the white noise.

  • keep_parity (bool, optional) – Whether keep parity of dimensions after padding.

  • full_output (bool, optional) – Whether to also return the indices of input frame within the padded frame (in addition to padded frame).

Returns:

  • array_out (numpy ndarray) – Padded array.

  • ori_indices (tuple) – [returned if full_output=True] Indices of the bottom left and top right vertices of the original image within the padded array: (y0, yN, x0, xN).

vip_hci.preproc.derotation module

Module with frame de-rotation routine for ADI.

[LAR97] (1,2)
Larkin et al. 1997
Fast Fourier method for the accurate rotation of sampled images
Optics Communications, Volume 154, Issue 1-3, pp. 99-106
vip_hci.preproc.derotation.cube_derotate(array, angle_list, imlib='vip-fft', interpolation='lanczos4', cxy=None, nproc=1, border_mode='constant', mask_val=nan, edge_blend=None, interp_zeros=False, ker=1)[source]

Rotate a cube (3d array or image sequence) providing a vector or corresponding angles.

Serves for rotating an ADI sequence to a common north given a vector with the corresponding parallactic angles for each frame.

Parameters:
  • array (numpy ndarray) – Input 3d array, cube.

  • angle_list (list) – Vector containing the parallactic angles.

  • imlib (str, optional) – See the documentation of the vip_hci.preproc.frame_rotate function.

  • interpolation (str, optional) – See the documentation of the vip_hci.preproc.frame_rotate function.

  • cxy (tuple of int, optional) – Coordinates X,Y of the point with respect to which the rotation will be performed. By default the rotation is done with respect to the center of the frames, as it is returned by the function vip_hci.var.frame_center.

  • nproc (int, optional) – Whether to rotate the frames in the sequence in a multi-processing fashion. Only useful if the cube is significantly large (frame size and number of frames).

  • border_mode (str, optional) – See the documentation of the vip_hci.preproc.frame_rotate function.

  • mask_val (flt, opt) – See the documentation of the vip_hci.preproc.frame_rotate function.

  • edge_blend (str, optional) – See the documentation of the vip_hci.preproc.frame_rotate function.

  • interp_zeros (str, optional) – See the documentation of the vip_hci.preproc.frame_rotate function.

  • ker (int, optional) – See the documentation of the vip_hci.preproc.frame_rotate function.

Returns:

array_der – Resulting cube with de-rotated frames.

Return type:

numpy ndarray

vip_hci.preproc.derotation.frame_rotate(array, angle, imlib='vip-fft', interpolation='lanczos4', cxy=None, border_mode='constant', mask_val=nan, edge_blend=None, interp_zeros=False, ker=1)[source]

Rotate a frame or 2D array.

Parameters:
  • array (numpy ndarray) – Input image, 2d array.

  • angle (float) – Rotation angle.

  • imlib ({'opencv', 'skimage', 'vip-fft'}, str optional) – Library used for image transformations. Opencv is faster than skimage or ‘vip-fft’, but vip-fft slightly better preserves the flux in the image (followed by skimage with a biquintic interpolation). ‘vip-fft’ corresponds to the FFT-based rotation method described in [LAR97], and implemented in this module. Best results are obtained with images without any sharp intensity change (i.e. no numerical mask). Edge-blending and/or zero-interpolation may help if sharp transitions are unavoidable.

  • interpolation (str, optional) – [Only used for imlib=’opencv’ or imlib=’skimage’] For Skimage the options are: ‘nearneig’, bilinear’, ‘biquadratic’, ‘bicubic’, ‘biquartic’ or ‘biquintic’. The ‘nearneig’ interpolation is the fastest and the ‘biquintic’ the slowest. The ‘nearneig’ is the poorer option for interpolation of noisy astronomical images. For Opencv the options are: ‘nearneig’, ‘bilinear’, ‘bicubic’ or ‘lanczos4’. The ‘nearneig’ interpolation is the fastest and the ‘lanczos4’ the slowest and more accurate. ‘lanczos4’ is the default for Opencv and ‘biquartic’ for Skimage.

  • cxy (float, optional) – Coordinates X,Y of the point with respect to which the rotation will be performed. By default the rotation is done with respect to the center of the frame.

  • border_mode ({'constant', 'edge', 'symmetric', 'reflect', 'wrap'}, str opt) – Pixel extrapolation method for handling the borders. ‘constant’ for padding with zeros. ‘edge’ for padding with the edge values of the image. ‘symmetric’ for padding with the reflection of the vector mirrored along the edge of the array. ‘reflect’ for padding with the reflection of the vector mirrored on the first and last values of the vector along each axis. ‘wrap’ for padding with the wrap of the vector along the axis (the first values are used to pad the end and the end values are used to pad the beginning). Default is ‘constant’.

  • mask_val (flt, opt) – If any numerical mask in the image to be rotated, what are its values? Will only be used if a strategy to mitigate Gibbs effects is adopted - see below.

  • edge_blend (str, opt {None,'noise','interp','noise+interp'}) –

    Whether to blend the edges, by padding nans then inter/extrapolate them with a gaussian filter. Slower but can significantly reduce ringing artefacts from Gibbs phenomenon, in particular if several consecutive rotations are involved in your image processing.

    • ’noise’: pad with small amplitude noise inferred from neighbours

    • ’interp’: interpolated from neighbouring pixels using Gaussian kernel.

    • ’noise+interp’: sum both components above at masked locations.

    Original mask will be placed back after rotation.

  • interp_zeros (bool, opt) – [only used if edge_blend is not None] Whether to interpolate zeros in the frame before (de)rotation. Not dealing with them can induce a Gibbs phenomenon near their location. However, this flag should be false if rotating a binary mask.

  • ker (float, opt) – Size of the Gaussian kernel used for interpolation.

Returns:

array_out – Resulting frame.

Return type:

numpy ndarray

vip_hci.preproc.derotation.rotate_fft(array, angle)[source]

Rotate a frame or 2D array using Fourier transforms.

Rotation is equivalent to 3 consecutive linear shears, or 3 consecutive 1D FFT phase shifts. See details in [LAR97].

Parameters:
  • array (numpy ndarray) – Input image, 2d array.

  • angle (float) – Rotation angle.

Returns:

array_out – Resulting frame.

Return type:

numpy ndarray

Note

This method is slower than interpolation methods (e.g. opencv/lanczos4 or ndimage), but preserves the flux better (by construction it preserves the total power). It is more prone to large-scale Gibbs artefacts, so make sure no sharp edge nor bad pixels are present in the image to be rotated.

Note

Warning: if input frame has even dimensions, the center of rotation will NOT be between the 4 central pixels, instead it will be on the top right of those 4 pixels. Make sure your images are centered with respect to that pixel before rotation.

vip_hci.preproc.parangles module

Module with frame parallactic angles calculations and de-rotation routines for ADI.

[MEE98]
J. Meeus 1998
Astronomical Algorithms
Book
vip_hci.preproc.parangles.check_pa_vector(angle_list, unit='deg')[source]

Check if the angle list has the right format to avoid any bug in the psfsub algorithms.

The right format complies to 3 criteria:

  1. angles are expressed in degree

  2. the angles are positive

  3. there is no jump of more than 180 deg between consecutive values (e.g. no jump like [350deg,355deg,0deg,5deg] => replaced by [350deg,355deg,360deg,365deg])

Parameters:
  • angle_list (1D-numpy ndarray) – Vector containing the derotation angles

  • unit (string, {'deg','rad'}, optional) – The unit type of the input angle list

Returns:

angle_list – Vector containing the derotation angles (after correction to comply with the 3 criteria, if needed)

Return type:

1-D numpy ndarray

vip_hci.preproc.parangles.compute_derot_angles_cd(objname_tmp_A, digit_format=3, objname_tmp_B='', inpath='./', skew=False, writing=False, outpath='./', list_obj=None, cd11_key='CD1_1', cd12_key='CD1_2', cd21_key='CD2_1', cd22_key='CD2_2', verbose=False)[source]

Function that returns a numpy vector of angles to derotate datacubes so as to match North up, East left, based on the CD matrix information contained in the header. In case the PosAng keyword is present in the header and there is no skewness between x and y axes, favor the use of function compute_derot_angles_PA (more precise as it averages for the middle of the exposure). The output is in appropriate format for the pca algorithm in the sense that:

  1. all angles of the output are in degrees

  2. all angles of the ouput are positive

  3. there is no jump of more than 180 deg between consecutive values (e.g. no jump like [350deg,355deg,0deg,5deg] => replaced by [350deg,355deg,360deg,365deg])

Parameters:
  • objname_tmp_A (string) – Contains the common name of the cubes BEFORE the digits

  • digit_format (int, optional) – Number of digits in the name of the cube. The digits are supposed to be the only changing part in the name of one cube to another.

  • objname_tmp_B (string, optional) – Contains the name of the cubes AFTER the digits

  • inpath (string, optional) – Contains the full path of the directory with the data

  • skew (bool, optional) – True if you know there is a different rotation between y- and x- axes. The code also detects automatically if there is >1deg skew between y and x axes. In case of skewing, 2 vectors of derotation angles are returned: one for x and one for y, instead of only one vector.

  • writing (bool, optional) – True if you want to write the derotation angles in a txt file.

  • outpath (string, opt) – Contains the full path of the directory where you want the txt file to be saved.

  • list_obj (integer list or 1-D array or None, optional) – List of the digits corresponding to the cubes to be considered. If not provided, the function will consider automatically all the cubes with objname_tmp_A+digit+objname_tmp_B+’.fits’ name structure in the provided “inpath”.

  • cd11_key (strings, optional) –

    Name of the keywords to be looked up in the header, to provide the:

    • partial of first axis coordinate w.r.t. x (cd11_key)

    • partial of first axis coordinate w.r.t. y (cd12_key)

    • partial of second axis coordinate w.r.t. x (cd21_key)

    • partial of second axis coordinate w.r.t. y (cd22_key)

    Default values are the ones in the headers of ESO or HST fits files. For more information, go to: http://www.stsci.edu/hst/HST_overview/documents/multidrizzle/ch44.html

  • cd12_key (strings, optional) –

    Name of the keywords to be looked up in the header, to provide the:

    • partial of first axis coordinate w.r.t. x (cd11_key)

    • partial of first axis coordinate w.r.t. y (cd12_key)

    • partial of second axis coordinate w.r.t. x (cd21_key)

    • partial of second axis coordinate w.r.t. y (cd22_key)

    Default values are the ones in the headers of ESO or HST fits files. For more information, go to: http://www.stsci.edu/hst/HST_overview/documents/multidrizzle/ch44.html

  • cd21_key (strings, optional) –

    Name of the keywords to be looked up in the header, to provide the:

    • partial of first axis coordinate w.r.t. x (cd11_key)

    • partial of first axis coordinate w.r.t. y (cd12_key)

    • partial of second axis coordinate w.r.t. x (cd21_key)

    • partial of second axis coordinate w.r.t. y (cd22_key)

    Default values are the ones in the headers of ESO or HST fits files. For more information, go to: http://www.stsci.edu/hst/HST_overview/documents/multidrizzle/ch44.html

  • cd22_key (strings, optional) –

    Name of the keywords to be looked up in the header, to provide the:

    • partial of first axis coordinate w.r.t. x (cd11_key)

    • partial of first axis coordinate w.r.t. y (cd12_key)

    • partial of second axis coordinate w.r.t. x (cd21_key)

    • partial of second axis coordinate w.r.t. y (cd22_key)

    Default values are the ones in the headers of ESO or HST fits files. For more information, go to: http://www.stsci.edu/hst/HST_overview/documents/multidrizzle/ch44.html

  • verbose (bool, optional) – True if you want more info to be printed.

Examples

If your cubes are: /home/foo/out_cube_obj_HK_025_000_sorted.fits, /home/foo/out_cube_obj_HK_025_001_sorted.fits, /home/foo/out_cube_obj_HK_025_002_sorted.fits, etc, the arguments should be:

objname_tmp_A = 'out_cube_obj_HK_025_'
digit_format = 3
objname_tmp_B = '_sorted'
inpath = '/home/foo/'
Returns:

angle_list – vector of angles corresponding to the angular difference between the positive y axis and the North in the image. sign convention: positive angles in anti-clockwise direction. Opposite values are applied when rotating the image to match North up. Note: if skew is set to True, there are 2 angle_list vectors returned; the first to rotate the x-axis and the second for the y-axis.

Return type:

1-D numpy ndarray

vip_hci.preproc.parangles.compute_derot_angles_pa(objname_tmp_A, digit_format=3, objname_tmp_B='', inpath='./', writing=False, outpath='./', list_obj=None, PosAng_st_key='HIERARCH ESO ADA POSANG', PosAng_nd_key='HIERARCH ESO ADA POSANG END', verbose=False)[source]

Function that returns a numpy vector of angles to derotate datacubes so as to match North up, East left, based on the mean of the Position Angle at the beginning and the end of the exposure. => It is twice more precise than function derot_angles_CD (there can be >1deg difference in the resulting angle vector returned for fast rotators with long exposures!), but it requires:

  1. a header keyword for both the position angle at start and end of exposure

  2. no skewness of the frames

The output is in appropriate format for the pca algorithm in the sense that:

  1. all angles of the output are in degrees

  2. all angles of the ouput are positive

  3. there is no jump of more than 180 deg between consecutive values (e.g. no jump like [350deg,355deg,0deg,5deg] => replaced by [350deg,355deg,360deg,365deg])

Parameters:
  • objname_tmp_A (string) – Contains the common name of the cubes BEFORE the digits

  • digit_format (int, optional) – Number of digits in the name of the cube. The digits are supposed to be the only changing part in the name of one cube to another.

  • objname_tmp_B (string, optional) – Contains the name of the cubes AFTER the digits

  • inpath (string, optional) – Contains the full path of the directory with the data

  • writing (bool, optional) – True if you want to write the derotation angles in a txt file.

  • outpath (string, optional) – Contains the full path of the directory where you want the txt file to be saved.

  • list_obj (integer list or 1-D array, optional) – List of the digits corresponding to the cubes to be considered. If not provided, the function will consider automatically all the cubes with objname_tmp_A+digit+objname_tmp_B+’.fits’ name structure in the provided “inpath”.

  • PosAng_st_key (strings, optional) – Name of the keywords to be looked up in the header, to provide the PA from North at start and end of integration.

  • PosAng_nd_key (strings, optional) – Name of the keywords to be looked up in the header, to provide the PA from North at start and end of integration.

  • verbose (bool, optional) – True if you want more info to be printed.

Examples

If your cubes are: /home/foo/out_cube_obj_HK_025_000_sorted.fits, /home/foo/out_cube_obj_HK_025_001_sorted.fits, /home/foo/out_cube_obj_HK_025_002_sorted.fits, etc, the arguments should be:

objname_tmp_A = 'out_cube_obj_HK_025_'
digit_format = 3
objname_tmp_B = '_sorted'
inpath = '/home/foo/'
Returns:

angle_list – vector of angles corresponding to the angular difference between the positive y axis and the North in the image. sign convention: positive angles in anti-clockwise direction. Opposite values are applied when rotating the image to match North up.

Return type:

1-D numpy ndarray

vip_hci.preproc.parangles.compute_paral_angles(header, latitude, ra_key, dec_key, lst_key, acqtime_key, date_key='DATE-OBS')[source]

Calculates the parallactic angle for a frame, taking coordinates and local sidereal time from fits-headers (frames taken in an alt-az telescope with the image rotator off).

The coordinates in the header are assumed to be J2000 FK5 coordinates. The spherical trigonometry formula for calculating the parallactic angle is taken from [MEE98].

Parameters:
  • header (dictionary) – Header of current frame.

  • latitude (float) – Latitude of the observatory in degrees. The dictionaries in vip_hci.conf.param can be used like: latitude=LBT[‘latitude’].

  • ra_key (strings) – Keywords where the values are stored in the header.

  • dec_key (strings) – Keywords where the values are stored in the header.

  • lst_key (strings) – Keywords where the values are stored in the header.

  • acqtime_key (strings) – Keywords where the values are stored in the header.

  • date_key (strings) – Keywords where the values are stored in the header.

Returns:

pa.value – Parallactic angle in degrees for current header (frame).

Return type:

float

vip_hci.preproc.recentering module

Module containing functions for cubes frame registration.

[GUI08] (1,2)
Guizar-Sicairos et al. 2008
Efficient subpixel image registration algorithms
Optics Letters, Volume 33, Issue 2, p. 156
[PUE15] (1,2)
Pueyo et al. 2015
Reconnaissance of the HR 8799 Exosolar System. II. Astrometry and Orbital Motion
The Astrophysical Journal, Volume 803, Issue 1, p. 31
vip_hci.preproc.recentering.cube_recenter_2dfit(array, xy=None, fwhm=4, subi_size=5, model='gauss', nproc=1, imlib='vip-fft', interpolation='lanczos4', offset=None, negative=False, threshold=False, sigfactor=2, fix_neg=False, params_2g=None, border_mode='reflect', save_shifts=False, full_output=False, verbose=True, debug=False, plot=True)[source]

Recenter the frames of a cube such that the centroid of a fitted 2D model is placed at the center of the frames.

The shifts are found by fitting a 2D Gaussian, Moffat, Airy or double Gaussian (positive+negative), as set by model, to a subimage centered at xy. This assumes the frames don’t have too large shifts (<5px). The frames are then shifted using the function frame_shift().

Parameters:
  • array (numpy ndarray) – Input cube.

  • xy (tuple of integers or floats) – Integer coordinates of the center of the subimage. For the double Gaussian fit with fixed negative gaussian (fix_neg=True), this should correspond to the exact location of the center of the negative Gaussian (e.g. the center of the coronagraph mask). In that case a tuple of floats is accepted.

  • fwhm (float or numpy ndarray) – FWHM size in pixels, either one value (float) that will be the same for the whole cube, or an array of floats with the same dimension as the 0th dim of array, containing the fwhm for each channel (e.g. in the case of an IFS cube, where the FWHM varies with wavelength).

  • subi_size (int, optional) – Size of the square subimage sides in pixels.

  • model (str, optional) – Sets the type of fit to be used. ‘gauss’ for a 2d Gaussian fit, ‘moff’ for a 2d Moffat fit, ‘airy’ for a 2d Airy disk fit, and ‘2gauss’ for a 2d double Gaussian (positive+negative) fit.

  • nproc (int or None, optional) – Number of processes (>1) for parallel computing. If 1 then it runs in serial. If None the number of processes will be set to (cpu_count()/2).

  • imlib (str, optional) – See the documentation of the vip_hci.preproc.frame_shift function.

  • interpolation (str, optional) – See the documentation of the vip_hci.preproc.frame_shift function.

  • offset (tuple of floats, optional) – If None the region of the frames used for the 2d Gaussian/Moffat fit is shifted to the center of the images (2d arrays). If a tuple is given it serves as the offset of the fitted area wrt the center of the 2d arrays.

  • negative (bool, optional) – If True a negative 2d Gaussian/Moffat fit is performed.

  • fix_neg (bool, optional) – In case of a double gaussian fit, whether to fix the parameters of the megative gaussian. If True, they should be provided in params_2g.

  • params_2g (None or dictionary, optional) –

    In case of a double gaussian fit, dictionary with either fixed or first guess parameters for the double gaussian. E.g.: params_2g = {‘fwhm_neg’: 3.5, ‘fwhm_pos’: (3.5,4.2), ‘theta_neg’: 48., ‘theta_pos’:145., ‘neg_amp’: 0.5}

    • fwhm_neg: float or tuple with fwhm of neg gaussian

    • fwhm_pos: can be a tuple for x and y axes of pos gaussian (replaces fwhm)

    • theta_neg: trigonometric angle of the x axis of the neg gaussian (deg)

    • theta_pos: trigonometric angle of the x axis of the pos gaussian (deg)

    • neg_amp: amplitude of the neg gaussian wrt the amp of the positive one

    Note: it is always recommended to provide theta_pos and theta_neg for a better fit.

  • threshold (bool, optional) – If True the background pixels (estimated using sigma clipped statistics) will be replaced by small random Gaussian noise (recommended for 2g).

  • sigfactor (float, optional) – If thresholding is performed, set the the threshold in terms of gaussian sigma in the subimage (will depend on your cropping size).

  • border_mode ({'reflect', 'nearest', 'constant', 'mirror', 'wrap'}) – Points outside the boundaries of the input are filled accordingly. With ‘reflect’, the input is extended by reflecting about the edge of the last pixel. With ‘nearest’, the input is extended by replicating the last pixel. With ‘constant’, the input is extended by filling all values beyond the edge with zeros. With ‘mirror’, the input is extended by reflecting about the center of the last pixel. With ‘wrap’, the input is extended by wrapping around to the opposite edge. Default is ‘reflect’.

  • save_shifts (bool, optional) – Whether to save the shifts to a file in disk.

  • full_output (bool, optional) – Whether to return 2 1d arrays of shifts along with the recentered cube or not.

  • verbose (bool, optional) – Whether to print to stdout the timing or not.

  • debug (bool, optional) – If True the details of the fitting are shown. Won’t work when the cube contains >20 frames (as it might produce an extremely long output).

  • plot (bool, optional) – If True, the shifts are plotted.

Returns:

  • array_rec (numpy ndarray) – The recentered cube.

  • y (numpy ndarray) – [full_output=True] 1d array with the shifts in y.

  • x (numpy ndarray) – [full_output=True] 1d array with the shifts in x.

vip_hci.preproc.recentering.cube_recenter_dft_upsampling(array, upsample_factor=100, subi_size=None, center_fr1=None, negative=False, fwhm=4, imlib='vip-fft', interpolation='lanczos4', mask=None, border_mode='reflect', log=False, collapse='median', full_output=False, verbose=True, nproc=None, save_shifts=False, debug=False, plot=True, **collapse_args)[source]

Recenter a cube of frames using the DFT upsampling method as proposed in [GUI08] and implemented in the phase_cross_correlation function from scikit-image.

The algorithm (DFT upsampling) obtains an initial estimate of the cross-correlation peak by an FFT and then refines the shift estimation by upsampling the DFT only in a small neighborhood of that estimate by means of a matrix-multiply DFT.

Optionally, after alignment of all images to the first one, a 2D Gaussian fit can be made to the mean image to recenter them based on the location of the Gaussian centroid. This second stage is performed if subi_size is not None.

Parameters:
  • array (numpy ndarray) – Input cube.

  • upsample_factor (int, optional) – Upsampling factor (default 100). Images will be registered to within 1/upsample_factor of a pixel. The larger the slower the algorithm.

  • subi_size (int or None, optional) – Size of the square subimage sides in pixels, used to find the centroid of the mean aligned cube image (i.e. after DFT-based registration). If subi_size is None then the first frame is assumed to be centered already.

  • (cy_1 (center_fr1 =) – [subi_size != None] Coordinates of the center of the subimage for fitting a 2d Gaussian. Since the first part of the function of the algorithm aligns all subsequent frames to the first one. This tuple should be the rough coordinates of the centroid in the first frame. If not provided, the function considers the center of the images.

  • cx_1) (Tuple, optional) – [subi_size != None] Coordinates of the center of the subimage for fitting a 2d Gaussian. Since the first part of the function of the algorithm aligns all subsequent frames to the first one. This tuple should be the rough coordinates of the centroid in the first frame. If not provided, the function considers the center of the images.

  • negative (bool, optional) – [subi_size != None] If True the final centroiding is done with a negative 2D Gaussian fit, instead of a positive one.

  • fwhm (float, optional) – [subi_size != None] First guess of the FWHM in pixels for the Gaussian fit.

  • nproc (int or None, optional) – Number of processes (>1) for parallel computing. If 1 then it runs in serial. If None the number of processes will be set to (cpu_count()/2).

  • imlib (str, optional) – See the documentation of the vip_hci.preproc.frame_shift function.

  • interpolation (str, optional) – See the documentation of the vip_hci.preproc.frame_shift function.

  • mask (2D np.ndarray, optional) – Binary mask indicating where the cross-correlation should be calculated in the images. If provided, should be the same size as array frames. [Note: requires skimage >= 0.18.0]

  • border_mode ({'reflect', 'nearest', 'constant', 'mirror', 'wrap'}) – Points outside the boundaries of the input are filled accordingly. With ‘reflect’, the input is extended by reflecting about the edge of the last pixel. With ‘nearest’, the input is extended by replicating the last pixel. With ‘constant’, the input is extended by filling all values beyond the edge with zeros. With ‘mirror’, the input is extended by reflecting about the center of the last pixel. With ‘wrap’, the input is extended by wrapping around to the opposite edge. Default is ‘reflect’.

  • log (bool, optional) – Whether to run the cross-correlation algorithm on images converted in log scale. This can be useful to leverage the whole extent of the PSF and be less dominated by the brightest central pixels.

  • collapse ({'median', 'mean', 'sum', 'max', 'trimmean', 'absmean', 'wmean'}) – Method used to collapse the aligned cube before 2D Gaussian fit. Should be an argument accepted by the vip_hci.preproc.cube_collapse function.

  • full_output (bool, optional) – Whether to return 2 1d arrays of shifts along with the recentered cube or not.

  • verbose (bool, optional) – Whether to print to stdout the timing or not.

  • save_shifts (bool, optional) – Whether to save the shifts to a file in disk.

  • debug (bool, optional) – Whether to print to stdout the shifts or not.

  • plot (bool, optional) – If True, the shifts are plotted.

  • collapse_args (opt) – Additional options passed to the vip_hci.preproc.cube_collapse function.

Returns:

  • array_recentered (numpy ndarray) – The recentered cube.

  • y (numpy ndarray) – [full_output=True] 1d array with the shifts in y.

  • x (numpy ndarray) – [full_output=True] 1d array with the shifts in x.

Note

This function uses the implementation from scikit-image of the algorithm described in [GUI08]. This algorithm registers two images (2-D rigid translation) within a fraction of a pixel specified by the user. Instead of computing a zero-padded FFT (fast Fourier transform), this code uses selective upsampling by a matrix-multiply DFT (discrete FT) to dramatically reduce computation time and memory without sacrificing accuracy. With this procedure all the image points are used to compute the upsampled cross-correlation in a very small neighborhood around its peak.

vip_hci.preproc.recentering.cube_recenter_radon(array, full_output=False, verbose=True, imlib='vip-fft', interpolation='lanczos4', border_mode='reflect', nproc=None, **kwargs)[source]

Recenter a cube using the Radon transform, as in [PUE15].

The function loops through its frames, relying on the frame_center_radon function for the recentering.

Parameters:
  • array (numpy ndarray) – Input 3d array or cube.

  • full_output ({False, True}, bool optional) – If True the recentered cube is returned along with the y and x shifts.

  • verbose ({True, False}, bool optional) – Whether to print timing and intermediate information to stdout.

  • imlib (str, optional) – See the documentation of the vip_hci.preproc.frame_shift function.

  • interpolation (str, optional) – See the documentation of the vip_hci.preproc.frame_shift function.

  • border_mode ({'reflect', 'nearest', 'constant', 'mirror', 'wrap'}) – Points outside the boundaries of the input are filled accordingly. With ‘reflect’, the input is extended by reflecting about the edge of the last pixel. With ‘nearest’, the input is extended by replicating the last pixel. With ‘constant’, the input is extended by filling all values beyond the edge with zeros. With ‘mirror’, the input is extended by reflecting about the center of the last pixel. With ‘wrap’, the input is extended by wrapping around to the opposite edge. Default is ‘reflect’.

  • nproc (int, optional) – Number of processes for parallel computing. If None the number of processes will be set to cpu_count()/2.

  • kwargs – Other optional parameters for vip_hci.preproc.frame_center_radon function, such as cropsize, hsize, step, satspots_cfg, mask_center, hpf, filter_fwhm, nproc or debug.

Returns:

  • array_rec (3d ndarray) – Recentered cube.

  • y, x (1d arrays of floats) – [full_output] Shifts in y and x.

  • dyx (1d array of floats) – [full_output] Array of uncertainty on center in pixels.

vip_hci.preproc.recentering.cube_recenter_satspots(array, xy, subi_size=19, sigfactor=6, plot=True, fit_type='moff', lbda=None, filter_freq=(0, 0), border_mode='constant', debug=False, verbose=True, full_output=False)[source]

Recenter an image cube based on satellite spots (more details in `.

The function relies on frame_center_satspots to align each image of the cube individually (details in vip_hci.preproc.frame_center_satspots). The function returns the recentered image cube, abd can also plot the histogram of the shifts, and calculate its statistics. The latter can help to assess the dispersion of the star center by using waffle/satellite spots (like those in VLT/SPHERE images) and evaluate the uncertainty of the position of the center.

Parameters:
  • array (numpy ndarray, 3d) – Input cube.

  • xy (tuple of 4 tuples of 2 elements) – Tuple with coordinates X,Y of the 4 satellite spots. When the spots are in an X configuration, the order is the following: top-left, top-right, bottom-left and bottom-right. When the spots are in an + (plus-like) configuration, the order is the following: top, right, left, bottom. If wavelength vector is not provided, assumes all sat spots of the cube are at a similar location. If wavelength is provided, only coordinates of the sat spots in the first channel should be provided. The boxes location in other channels will be scaled accordingly.

  • subi_size (int, optional) – Size of subimage where the fitting is done.

  • sigfactor (int, optional) – The background pixels will be thresholded before fitting a 2d Gaussian to the data using sigma clipped statistics. All values smaller than (MEDIAN + sigfactor*STDDEV) will be replaced by small random Gaussian noise.

  • plot (bool, optional) – Whether to plot the shifts.

  • fit_type (str, optional {'gaus','moff'}) – Type of 2d fit to infer the centroid of the satellite spots.

  • lbda (1d array or list, opt) – Wavelength vector. If provided, the subimages will be scaled accordingly to follow the motion of the satellite spots.

  • filter_freq (tuple of 2 floats, optional) – If the first (resp. second) element of the tuple is larger than 0, a high-pass (resp. low-pass) filter is applied to the image, before fitting the satellite spots. The elements should correspond to the fwhm_size of the frame_filter_highpass and frame_filter_lowpass functions, respectively. If both elements are non-zero, both high-pass and low-pass filter of the image are applied, in that order. This can be useful to isolate the signal from the satellite spots.

  • border_mode ({'reflect', 'nearest', 'constant', 'mirror', 'wrap'}) – Points outside the boundaries of the input are filled accordingly. With ‘reflect’, the input is extended by reflecting about the edge of the last pixel. With ‘nearest’, the input is extended by replicating the last pixel. With ‘constant’, the input is extended by filling all values beyond the edge with zeros. With ‘mirror’, the input is extended by reflecting about the center of the last pixel. With ‘wrap’, the input is extended by wrapping around to the opposite edge. Default is ‘reflect’.

  • debug (bool, optional) – If True debug information is printed and plotted (fit and residuals, intersections and shifts). This has to be used carefully as it can produce too much output and plots.

  • verbose (bool, optional) – Whether to print to stdout the timing and additional info.

  • full_output (bool, optional) – Whether to return 2 1d arrays of shifts along with the recentered cube or not.

Returns:

  • array_rec – The shifted cube.

  • shift_y, shift_x – [full_output==True] Shifts Y,X to get to the true center for each image.

  • sat_y, sat_x – [full_output==True] Y,X positions of the satellite spots in each image. Order: top-left, top-right, bottom-left and bottom-right.

vip_hci.preproc.recentering.cube_recenter_via_speckles(cube_sci, cube_ref=None, alignment_iter=5, gammaval=1, min_spat_freq=0.5, max_spat_freq=3, fwhm=4, upsample_factor=100, debug=False, recenter_median=False, fit_type='gaus', negative=True, crop=True, subframesize=25, mask=None, ann_rad=0.5, ann_rad_search=False, ann_width=0.5, collapse='median', imlib='vip-fft', interpolation='lanczos4', border_mode='reflect', log=True, plot=True, full_output=False, nproc=None, **collapse_args)[source]

Register frames based on the median speckle pattern.

The function also optionally centers images based on the position of the vortex null in the median frame (through a negative Gaussian fit or an annulus fit). By default, images are filtered to isolate speckle spatial frequencies, and converted to log-scale before the cross-correlation based alignment. The image cube should already be centered within ~2px accuracy before being passed to this function (e.g. through an eyeball crop using vip_hci.preproc.cube_crop_frames).

Parameters:
  • cube_sci (numpy ndarray) – Science cube.

  • cube_ref (numpy ndarray) – Reference cube (e.g. for NIRC2 data in RDI mode).

  • alignment_iter (int, optional) – Number of alignment iterations (recomputes median after each iteration).

  • gammaval (int, optional) – Applies a gamma correction to emphasize speckles (useful for faint stars).

  • min_spat_freq (float, optional) – Spatial frequency for low pass filter.

  • max_spat_freq (float, optional) – Spatial frequency for high pass filter.

  • fwhm (float, optional) – Full width at half maximum.

  • upsample_factor (int, optional) – Upsampling factor (default 100). Images will be registered to within 1/upsample_factor of a pixel. The larger the slower the algorithm.

  • debug (bool, optional) – Outputs extra info.

  • recenter_median (bool, optional) – Recenter the frames at each iteration based on a 2d fit.

  • fit_type (str, optional) – If recenter_median is True, this is the model to which the image is fitted to for recentering. ‘gaus’ works well for NIRC2_AGPM data. ‘ann’ works better for NACO+AGPM data.

  • negative (bool, optional) – If True, uses a negative gaussian fit to determine the center of the median frame.

  • crop (bool, optional) – Whether to calculate the recentering on a cropped version of the cube that is speckle-dominated (recommended).

  • subframesize (int, optional) – Sub-frame window size used. Should cover the region where speckles are the dominant noise source.

  • mask (2D np.ndarray, optional) – Binary mask indicating where the cross-correlation should be calculated in the images. If provided, should be the same size as array frames.

  • ann_rad (float, optional) – [if fit_type=’ann’] The expected inner radius of the annulus in FWHM.

  • ann_rad_search (bool) – [if fit_type=’ann’] Whether to also search for optimal radius.

  • ann_width (float, optional) – [if fit_type=’ann’] The expected radial width of the annulus in FWHM.

  • collapse ({'median', 'mean', 'sum', 'max', 'trimmean', 'absmean', 'wmean'}) – Method used to collapse the aligned cube before 2D Gaussian fit. Should be an argument accepted by the vip_hci.preproc.cube_collapse function.

  • imlib (str, optional) – Image processing library to use.

  • interpolation (str, optional) – Interpolation method to use.

  • border_mode ({'reflect', 'nearest', 'constant', 'mirror', 'wrap'}) – Points outside the boundaries of the input are filled accordingly. With ‘reflect’, the input is extended by reflecting about the edge of the last pixel. With ‘nearest’, the input is extended by replicating the last pixel. With ‘constant’, the input is extended by filling all values beyond the edge with zeros. With ‘mirror’, the input is extended by reflecting about the center of the last pixel. With ‘wrap’, the input is extended by wrapping around to the opposite edge. Default is ‘reflect’.

  • log (bool) – Whether to run the cross-correlation algorithm on images converted in log scale. This can be useful to leverage the whole extent of the PSF and be less dominated by the brightest central pixels.

  • plot (bool, optional) – If True, the shifts are plotted.

  • full_output (bool, optional) – Whether to return more variables, useful for debugging.

  • **collapse_args – Additional arguments passed to the vip_hci.preproc.cube_collapse function.

Returns:

  • cube_reg_sci (numpy 3d ndarray) – Registered science cube

  • cube_reg_ref (numpy 3d ndarray) – [cube_ref!=None] Cube registered to science frames

  • cube_sci_lpf (numpy 3d ndarray) – [full_output=True] Low+high-pass filtered science cube

  • cube_stret (numpy 3d ndarray) – [full_output=True] cube_stret with stretched values used for cross-corr

  • cum_x_shifts_sci (numpy 1d array) – [full_output=True] Vector of x shifts for science frames

  • cum_y_shifts_sci (numpy 1d array) – [full_output=True] Vector of y shifts for science frames

  • cum_x_shifts_ref (numpy 1d array) – [full_output=True & cube_ref!=None] Vector of x shifts for ref frames

  • cum_y_shifts_ref (numpy 1d array) – [full_output=True & cube_ref!=None] Vector of y shifts for ref frames

vip_hci.preproc.recentering.cube_shift(cube, shift_y, shift_x, imlib='vip-fft', interpolation='lanczos4', border_mode='reflect', nproc=None)[source]

Shift the X-Y coordinates of a cube or 3D array by x and y values.

Parameters:
  • cube (numpy ndarray, 3d) – Input cube.

  • shift_y (float, list of floats or np.ndarray of floats) – Shifts in y and x directions for each frame. If the a single value is given then all the frames will be shifted by the same amount.

  • shift_x (float, list of floats or np.ndarray of floats) – Shifts in y and x directions for each frame. If the a single value is given then all the frames will be shifted by the same amount.

  • imlib (str, optional) – See the documentation of the vip_hci.preproc.frame_shift function.

  • interpolation (str, optional) – See the documentation of the vip_hci.preproc.frame_shift function.

  • border_mode (str, optional) – See the documentation of the vip_hci.preproc.frame_shift function.

  • nproc (int or None, optional) – Number of CPUs to use for multiprocessing. If None, will be automatically set to half the number of available CPUs.

Returns:

cube_out – Cube with shifted frames.

Return type:

numpy ndarray, 3d

vip_hci.preproc.recentering.frame_center_radon(array, cropsize=None, hsize_ini=1.0, step_ini=0.1, n_iter=5, tol=0.1, mask_center=None, nproc=None, satspots_cfg=None, theta_0=0, delta_theta=5, gauss_fit=True, hpf=True, filter_fwhm=8, imlib='vip-fft', interpolation='lanczos4', full_output=False, verbose=True, plot=True, debug=False)[source]

Find the center of a broadband (co-added) frame with speckles and satellite spots elongated towards the star (center).

The function uses the Radon transform implementation from scikit-image, and follow the algorithm presented in [PUE15].

Parameters:
  • array (numpy ndarray) – Input 2d array or image.

  • cropsize (None or odd int, optional) – Size in pixels of the cropped central area of the input array that will be used. It should be large enough to contain the bright elongated speckle or satellite spots.

  • hsize_ini (float, optional) – Size of the box for the grid search for first centering iteration. The frame is shifted to each direction from the center in a hsize length with a given step.

  • step_ini (float, optional) – The step of the coordinates change in the first step. Note: should not be too fine for efficiency as it is automatically refined at each step.

  • n_iter (int, optional) – Number of iterations for finer recentering. At each step, a finer step is considered based on the amplitude of the shifts found in the previous step. Iterations are particularly relevant when mask_center is not None, as the masked area will change from one iteration to the next.

  • tol (float, optional) – Absolute tolerance on relative shift from one iteration to the next to consider convergence. If the absolute value of the shift is found to be less than tol, the iterative algorithm is stopped.

  • mask_center (None or int, optional) – If None the central area of the frame is kept. If int a centered zero mask will be applied to the frame. By default the center isn’t masked.

  • nproc (int, optional) – Number of processes for parallel computing. If None the number of processes will be set to cpu_count()/2.

  • satspots_cfg (None or str ('x', '+' or 'custom'), opt) – If satellite spots are present, provide a string corresponding to the configuration of the satellite spots: as a cross (‘x’), as a plus sign (‘+’) or ‘custom’ (provide theta_0). Leave to None if no satellite spots present. Note: setting satspots_cfg to non-None value leads to varying performance depending on dataset.

  • theta_0 (float between [0,90[, optional) – Azimuth of the first satellite spot. Only considered if satspots_cfg is set to ‘custom’.

  • delta_theta (float, optional) – Azimuthal half-width in degrees of the slices considered along a ‘+’ or ‘x’ pattern to calculate the Radon transform. E.g. if set to 5 for ‘x’ configuration, it will consider slices from 40 to 50 deg in each quadrant.

  • hpf (bool, optional) – Whether to high-pass filter the images

  • filter_fwhm (float, optional) – In case of high-pass filtering, this is the FWHM of the low-pass filter used for subtraction to the original image to get the high-pass filtered image (i.e. should be >~ 2 x FWHM).

  • imlib (str, optional) – See the documentation of the vip_hci.preproc.frame_shift function.

  • interpolation (str, optional) – See the documentation of the vip_hci.preproc.frame_shift function.

  • full_output (bool, optional) – Whether to also return the cost map, and uncertainty on centering.

  • verbose (bool optional) – Whether to print to stdout some messages and info.

  • plot (bool, optional) – Whether to plot the radon cost function.

  • debug (bool, optional) – Whether to print and plot intermediate info.

Returns:

  • optimy, optimx (floats) – Values of the Y, X coordinates of the center of the frame based on the radon optimization. (always returned)

  • dxy (float) – [full_output=True] Uncertainty on center in pixels.

  • cost_bound (2d numpy array) – [full_output=True] Radon cost function surface.

vip_hci.preproc.recentering.frame_center_satspots(array, xy, subi_size=19, sigfactor=6, shift=False, imlib='vip-fft', interpolation='lanczos4', fit_type='moff', filter_freq=(0, 0), border_mode='reflect', debug=False, verbose=True)[source]

Find the center of a frame with satellite spots (relevant e.g. for VLT/SPHERE data).

The method used to determine the center is by centroiding the 4 spots via a 2D Gaussian fit and finding the intersection of the lines they create (see Notes). This method is very sensitive to the SNR of the satellite spots, therefore thresholding of the background pixels is performed. If the results are too extreme, the debug parameter will allow to see in depth what is going on with the fit (you may have to adjust the sigfactor for the background pixels thresholding).

Parameters:
  • array (numpy ndarray, 2d) – Image or frame.

  • xy (tuple of 4 tuples of 2 elements) – Tuple with coordinates X,Y of the 4 satellite spots. When the spots are in an X configuration, the order is the following: top-left, top-right, bottom-left and bottom-right. When the spots are in an + (cross-like) configuration, the order is the following: top, right, left, bottom.

  • subi_size (int, optional) – Size of subimage where the fitting is done.

  • sigfactor (int, optional) – The background pixels will be thresholded before fitting a 2d Gaussian to the data using sigma clipped statistics. All values smaller than (MEDIAN + sigfactor*STDDEV) will be replaced by small random Gaussian noise.

  • shift (bool, optional) – If True the image is shifted.

  • imlib (str, optional) – See the documentation of the vip_hci.preproc.frame_shift function.

  • interpolation (str, optional) – See the documentation of the vip_hci.preproc.frame_shift function.

  • fit_type (str, optional {'gaus','moff'}) – Type of 2d fit to infer the centroid of the satellite spots.

  • filter_freq (tuple of 2 floats, optional) – If the first (resp. second) element of the tuple is larger than 0, a high-pass (resp. low-pass) filter is applied to the image, before fitting the satellite spots. The elements should correspond to the fwhm_size of the frame_filter_highpass and frame_filter_lowpass functions, respectively. If both elements are non-zero, both high-pass and low-pass filter of the image are applied, in that order. This can be useful to better isolate the signal from the satellite spots.

  • border_mode ({'reflect', 'nearest', 'constant', 'mirror', 'wrap'}) – Points outside the boundaries of the input are filled accordingly. With ‘reflect’, the input is extended by reflecting about the edge of the last pixel. With ‘nearest’, the input is extended by replicating the last pixel. With ‘constant’, the input is extended by filling all values beyond the edge with zeros. With ‘mirror’, the input is extended by reflecting about the center of the last pixel. With ‘wrap’, the input is extended by wrapping around to the opposite edge. Default is ‘reflect’.

  • debug (bool, optional) – If True debug information is printed and plotted.

  • verbose (bool, optional) – If True the intersection and shifts information is printed out.

Returns:

  • array_rec (2d numpy array) – Shifted images. Only returned if ``shift=True``.

  • shifty, shiftx (floats) – Shift Y,X to get to the true center.

  • ceny, cenx (floats) – Center Y,X coordinates of the true center. Only returned if ``shift=True``.

Note

We are solving a linear system:

A1 * x + B1 * y = C1
A2 * x + B2 * y = C2

Cramer’s rule - solution can be found in determinants:

x = Dx/D
y = Dy/D

where D is main determinant of the system:

A1 B1
A2 B2

and Dx and Dy can be found from matrices:

C1 B1
C2 B2

and

A1 C1
A2 C2

C column consequently substitutes the coef. columns of x and y

L stores our coefs A, B, C of the line equations.

For D: L1[0] L1[1]   for Dx: L1[2] L1[1]   for Dy: L1[0] L1[2]
       L2[0] L2[1]           L2[2] L2[1]           L2[0] L2[2]
vip_hci.preproc.recentering.frame_shift(array, shift_y, shift_x, imlib='vip-fft', interpolation='lanczos4', border_mode='reflect')[source]

Shift a 2D array by shift_y, shift_x.

Parameters:
  • array (numpy ndarray) – Input 2d array.

  • shift_y (float) – Shifts in y and x directions.

  • shift_x (float) – Shifts in y and x directions.

  • imlib ({'opencv', 'ndimage-fourier', 'ndimage-interp', 'vip-fft'}, str opt) – Library or method used for performing the image shift. ‘ndimage-fourier’ or ‘vip-fft’: does a fourier shift operation and preserves better the pixel values - therefore the flux and photometry (wrapper of scipy.ndimage.fourier_shift). Interpolation-based shift (‘opencv’ and ‘ndimage-interp’) is faster but less accurate than the fourier shift. ‘opencv’ is recommended when speed is critical.

  • interpolation (str, optional) – Only used in case of imlib is set to ‘opencv’ or ‘ndimage-interp’ (Scipy.ndimage), where the images are shifted via interpolation. For Scipy.ndimage the options are: ‘nearneig’, bilinear’, ‘biquadratic’, ‘bicubic’, ‘biquartic’ or ‘biquintic’. The ‘nearneig’ interpolation is the fastest and the ‘biquintic’ the slowest. The ‘nearneig’ is the poorer option for interpolation of noisy astronomical images. For Opencv the options are: ‘nearneig’, ‘bilinear’, ‘bicubic’ or ‘lanczos4’. The ‘nearneig’ interpolation is the fastest and the ‘lanczos4’ the slowest and accurate. ‘lanczos4’ is the default for Opencv and ‘biquartic’ for Scipy.ndimage.

  • border_mode ({'reflect', 'nearest', 'constant', 'mirror', 'wrap'}) – For ‘opencv’ and ‘ndimage-interp’, points outside the boundaries of the input are filled according to the value of this parameter. With ‘reflect’, the input is extended by reflecting about the edge of the last pixel. With ‘nearest’, the input is extended by replicating the last pixel. With ‘constant’, the input is extended by filling all values beyond the edge with zeros. With ‘mirror’, the input is extended by reflecting about the center of the last pixel. With ‘wrap’, the input is extended by wrapping around to the opposite edge. Default is ‘reflect’. Note: for ‘ndimage-fourier’ default is ‘wrap’ (impossible to change), while border_mode is ‘constant’ (zeros) for ‘vip-fft’.

Returns:

array_shifted – Shifted 2d array.

Return type:

numpy ndarray

vip_hci.preproc.rescaling module

Module with frame px resampling/rescaling functions.

vip_hci.preproc.rescaling.check_scal_vector(scal_vec)[source]

Checks that the scaling factor list has the right format (i.e. all factors >= 1). If not, it returns the vector after normalization by the minimum value.

Parameters:

scal_vec (1d array or list) – Vector with the wavelengths.

Returns:

scal_vec – Vector containing the scaling factors (after correction to comply with the condition >= 1).

Return type:

numpy ndarray, 1d

vip_hci.preproc.rescaling.cube_px_resampling(array, scale, imlib='vip-fft', interpolation='lanczos4', keep_center=False, verbose=True)[source]

Resample the frames of a cube with a single scale factor. Can deal with NaN values.

Wrapper of frame_px_resampling. Useful when we need to upsample (upscaling) or downsample (pixel binning) a set of frames, e.g. an ADI cube.

Parameters:
  • array (3d numpy ndarray) – Input cube, 3d array.

  • scale (int, float or tuple) – Scale factor for upsampling or downsampling the frames in the cube. If a tuple it corresponds to the scale along x and y.

  • imlib (str, optional) – See the documentation of the vip_hci.preproc.frame_px_resampling function.

  • interpolation (str, optional) – See the documentation of the vip_hci.preproc.frame_px_resampling function.

  • keep_center (bool, opt) – If input dimensions are even and the star centered (i.e. on dim//2, dim//2), whether to keep the star centered after scaling, i.e. on (new_dim//2, new_dim//2). For a non-centered input cube, better to leave it to False.

  • verbose (bool, optional) – Whether to print out additional info such as the new cube shape.

Returns:

array_resc – Output cube with resampled frames.

Return type:

numpy ndarray

vip_hci.preproc.rescaling.cube_rescaling(array, scaling_list, ref_xy=None, imlib='vip-fft', interpolation='lanczos4', scaling_y=None, scaling_x=None)[source]

Rescale a cube by factors from scaling_list wrt a position.

Parameters:
  • array (numpy ndarray) – Input 3d array, cube.

  • scaling_list (1D numpy ndarray or list) – Scale corresponding to each frame in the cube.

  • ref_xy (float, optional) – Coordinates X,Y of the point with respect to which the rescaling will be performed. By default the rescaling is done with respect to the center of the frames; central pixel if the frames have odd size.

  • imlib (str optional) – See the documentation of vip_hci.preproc.cube_rescaling_wavelengths.

  • interpolation (str, optional) – See the documentation of vip_hci.preproc.cube_rescaling_wavelengths.

  • scaling_y (1D-array or list) – Scaling factor only for y axis. If provided, it takes priority on scaling_list.

  • scaling_x (1D-array or list) – Scaling factor only for x axis. If provided, it takes priority on scaling_list.

Returns:

array_sc – Resulting cube with rescaled frames.

Return type:

numpy ndarray

vip_hci.preproc.rescaling.cube_rescaling_wavelengths(cube, scal_list, full_output=True, inverse=False, y_in=None, x_in=None, imlib='vip-fft', interpolation='lanczos4', collapse='median', pad_mode='reflect')[source]

Scale/Descale a cube by scal_list, with padding. Can deal with NaN values.

Wrapper to scale or descale a cube by factors given in scal_list, without any loss of information (zero-padding if scaling > 1). Important: in case of IFS data, the scaling factors in scal_list should be >= 1 (ie. provide the scaling factors as for scaling to the longest wavelength channel).

Parameters:
  • cube (3D-array) – Data cube with frames to be rescaled.

  • scal_list (1D-array) – Vector of same dimension as the first dimension of datacube, containing the scaling factor for each frame.

  • full_output (bool, optional) – Whether to output just the rescaled cube (False) or also its median, the new y and x shapes of the cube, and the new centers cy and cx of the frames (True).

  • inverse (bool, optional) – Whether to inverse the scaling factors in scal_list before applying them or not; i.e. True is to descale the cube (typically after a first scaling has already been done)

  • y_in (int) – Initial y and x sizes, required for inverse=True. In case the cube is descaled, these values will be used to crop back the cubes/frames to their original size.

  • x_in (int) – Initial y and x sizes, required for inverse=True. In case the cube is descaled, these values will be used to crop back the cubes/frames to their original size.

  • imlib ({'opencv', 'ndimage', 'vip-fft'}, str optional) – Library used for image transformations. Opencv is faster than ndimage or skimage. ‘vip-fft’ corresponds to a FFT-based rescaling algorithm implemented in VIP (vip_hci.preproc.scale_fft).

  • interpolation (str, optional) – For ‘ndimage’ library: ‘nearneig’, bilinear’, ‘bicuadratic’, ‘bicubic’, ‘biquartic’, ‘biquintic’. The ‘nearneig’ interpolation is the fastest and the ‘biquintic’ the slowest. The ‘nearneig’ is the poorer option for interpolation of noisy astronomical images. For ‘opencv’ library: ‘nearneig’, ‘bilinear’, ‘bicubic’, ‘lanczos4’. The ‘nearneig’ interpolation is the fastest and the ‘lanczos4’ the slowest and accurate. ‘lanczos4’ is the default.

  • collapse ({'median', 'mean', 'sum', 'trimmean'}, str optional) – Sets the way of collapsing the frames for producing a final image.

  • pad_mode (str, optional) –

    One of the following string values:

    'constant'

    pads with a constant value

    'edge'

    pads with the edge values of array

    'linear_ramp'

    pads with the linear ramp between end_value and the array edge value.

    'maximum'

    pads with the maximum value of all or part of the vector along each axis

    'mean'

    pads with the mean value of all or part of the vector along each axis

    'median'

    pads with the median value of all or part of the vector along each axis

    'minimum'

    pads with the minimum value of all or part of the vector along each axis

    'reflect'

    pads with the reflection of the vector mirrored on the first and last values of the vector along each axis

    'symmetric'

    pads with the reflection of the vector mirrored along the edge of the array

    'wrap'

    pads with the wrap of the vector along the axis. The first values are used to pad the end and the end values are used to pad the beginning

Returns:

  • frame (2d array) – [full_output=False] The median of the rescaled cube.

  • cube (3d array) – [full_output=True] Rescaled cube

  • frame (2d array) – [full_output=True] Median of the rescaled cube – note it is in 2nd position if full_output is set to True.

  • y,x,cy,cx (floats) – [full_output=True] New y and x shapes of the cube, and the new centers cy and cx of the frames

vip_hci.preproc.rescaling.find_scal_vector(cube, lbdas, fluxes, mask=None, nfp=2, fm='stddev', simplex_options=None, debug=False, imlib='vip-fft', interpolation='lanczos4', hpf=False, fwhm_max=5, **kwargs)[source]

Find the optimal scaling factor for the channels of an IFS cube (or of dual-band pairs of images).

The algorithm finds the optimal scaling factor that minimizes residuals in the rescaled frames. It takes the inverse of the wavelength vector as a first guess, and uses a similar method as the negative fake companion technique, but minimizing residuals in either a mask or the whole field.

Parameters:
  • cube (3D-array) – Data cube with frames to be rescaled.

  • lbdas (1d array) – Vector with the wavelengths, used for first guess on scaling factor.

  • fluxes (1d array) – Vector with the (unsaturated) fluxes at the different wavelengths, used for first guess on flux factor.

  • mask (2D-array, opt) – Binary mask, with ones where the residual intensities should be evaluated. If None is provided, the whole field is used.

  • nfp (int, opt, {1,2}) – Number of free parameters: spatial scaling alone or spatial scaling + flux scaling.

  • fm (str, opt, {"sum","stddev"}) – Figure of merit to use: sum of squared residuals or stddev of residual pixels.

  • options (dict, optional) – The scipy.optimize.minimize options.

  • hpf (bool, optional) – Whether to high-pass filter images before searching for optimal scaling factors. Can help to not be biased by a diffuse halo, and just leverage speckle expansion.

  • fwhm_max (float, optional) – Maximum FWHM of the PSF across all wavelengths, in pixels. Only used if hpf is set to True.

  • **kwargs (optional) – Optional arguments to the scipy.optimize.minimize function

Returns:

  • scal_vec (numpy ndarray, 1d) – Vector containing the scaling factors (after correction to comply with the condition >= 1).

  • flux_vec (numpy ndarray, 1d [only returned if nfp==2]) – Vector containing the associated flux factors.

vip_hci.preproc.rescaling.frame_px_resampling(array, scale, imlib='vip-fft', interpolation='lanczos4', keep_center=False, verbose=False)[source]

Resample the pixels of a frame changing the frame size. Can deal with NaN values.

If scale < 1 then the frame is downsampled and if scale > 1 then its pixels are upsampled.

Warning: if imlib is not ‘vip-fft’, the input size is even and keep_center set to True, an additional interpolation (shifting by (0.5,0.5)px) may occur after rescaling, to ensure center location stays on (dim//2,dim//2).

Parameters:
  • array (numpy ndarray) – Input frame, 2d array.

  • scale (int, float or tuple) – Scale factor for upsampling or downsampling the frame. If a tuple it corresponds to the scale along x and y.

  • imlib ({'ndimage', 'opencv', 'vip-fft'}, optional) – Library used for image transformations. ‘vip-fft’ corresponds to a FFT-based rescaling algorithm implemented in VIP (vip_hci.preproc.scale_fft).

  • interpolation (str, optional) – For ‘ndimage’ library: ‘nearneig’, bilinear’, ‘biquadratic’, ‘bicubic’, ‘biquartic’, ‘biquintic’. The ‘nearneig’ interpolation is the fastest and the ‘biquintic’ the slowest. The ‘nearneig’ is the worst option for interpolation of noisy astronomical images. For ‘opencv’ library: ‘nearneig’, ‘bilinear’, ‘bicubic’, ‘lanczos4’. The ‘nearneig’ interpolation is the fastest and the ‘lanczos4’ the slowest and accurate.

  • keep_center (bool, opt) – If input dimensions are even and the star centered (i.e. on dim//2, dim//2), whether to keep the star centered after scaling, i.e. on (new_dim//2, new_dim//2). For a non-centered input frame, better to leave it to False.

  • verbose (bool, optional) – Whether to print out additional info such as the new image shape.

Returns:

array_resc – Output resampled frame.

Return type:

numpy ndarray

vip_hci.preproc.rescaling.frame_rescaling(array, ref_xy=None, scale=1.0, imlib='vip-fft', interpolation='lanczos4', scale_y=None, scale_x=None)[source]

Rescale a frame by a factor wrt a reference point.

The reference point is by default the center of the frame (typically the exact location of the star). However, it keeps the same dimensions.

Parameters:
  • array (numpy ndarray) – Input frame, 2d array.

  • ref_xy (float, optional) – Coordinates X,Y of the point wrt which the rescaling will be applied. By default the rescaling is done with respect to the center of the frame.

  • scale (float) – Scaling factor. If > 1, it will upsample the input array equally along y and x by this factor.

  • imlib ({'ndimage', 'opencv', 'vip-fft'}, optional) – Library used for image transformations. ‘vip-fft’ corresponds to a FFT-based rescaling algorithm implemented in VIP (vip_hci.preproc.scale_fft).

  • interpolation (str, optional) – For ‘ndimage’ library: ‘nearneig’, bilinear’, ‘biquadratic’, ‘bicubic’, ‘biquartic’, ‘biquintic’. The ‘nearneig’ interpolation is the fastest and the ‘biquintic’ the slowest. The ‘nearneig’ is the worst option for interpolation of noisy astronomical images. For ‘opencv’ library: ‘nearneig’, ‘bilinear’, ‘bicubic’, ‘lanczos4’. The ‘nearneig’ interpolation is the fastest and the ‘lanczos4’ the slowest and accurate.

  • scale_y (float) – Scaling factor only for y axis. If provided, it takes priority on scale parameter.

  • scale_x (float) – Scaling factor only for x axis. If provided, it takes priority on scale parameter.

Returns:

array_out – Resulting frame.

Return type:

numpy ndarray

vip_hci.preproc.rescaling.scale_fft(array, scale, ori_dim=False)[source]

Resample the frames of a cube with a single scale factor using a FFT-based method.

Parameters:
  • array (3d numpy ndarray) – Input cube, 3d array.

  • scale (int or float) – Scale factor for upsampling or downsampling the frames in the cube. If a tuple it corresponds to the scale along x and y.

  • ori_dim (bool, opt) – Whether to crop/pad scaled array in order to have the output with the same dimensions as the input array. By default, the x,y dimensions of the output are the closest integer to scale*dim_input, with the same parity as the input.

Returns:

array_resc – Output cube with resampled frames.

Return type:

numpy ndarray

vip_hci.preproc.skysubtraction module

Module with sky subtraction function.

[GOM17]
Gomez-Gonzalez et al. 2017
VIP: Vortex Image Processing Package for High-contrast Direct Imaging
The Astronomical Journal, Volume 154, p. 7
[HUN18]
Hunziker et al. 2018
PCA-based approach for subtracting thermal background emission in high-contrast imaging data
Astronomy & Astrophysics, Volume 611, p. 23
[REN23]
Ren 2023
Karhunen-Loève data imputation in high-contrast imaging
Astronomy & Astrophysics, Volume 679, p. 8
vip_hci.preproc.skysubtraction.cube_subtract_sky_pca(sci_cube, sky_cube, masks, ref_cube=None, ncomp=2, full_output=False)[source]

PCA-based sky subtraction as explained in [REN23] (see also [GOM17] and [HUN18]).

Parameters:
  • sci_cube (numpy ndarray) – 3d array of science frames.

  • sky_cube (numpy ndarray) – 3d array of sky frames.

  • masks (tuple of two numpy ndarray or one signle numpy ndarray) – Mask indicating the boat and anchor regions for the analysis. If two masks are provided, they will be assigned to mask_anchor and mask_boat in that order. If only one mask is provided, it will be used as the anchor, and the boat images will not be masked (i.e., full frames used).

  • ref_cube (numpy ndarray or None, opt) – Reference cube.

  • ncomp (int, opt) – Sets the number of PCs you want to use in the sky subtraction.

  • full_output (bool, opt) – Whether to also output pcs, reconstructed cube, residuals cube and derotated residual cube.

Notes

Masks can be created with the function vip_hci.var.create_ringed_spider_mask or vip_hci.var.get_annulus_segments (see Usage Example below).

Returns:

  • sci_cube_skysub (numpy ndarray) – Sky-subtracted science cube

  • ref_cube_skysub (numpy ndarray) – [If ref_cube is not None] Also returns sky-subtracted reference cube.

  • If full_output is set to True, returns (in the following order)

    • sky-subtracted science cube,

    • sky-subtracted reference cube (if any provided),

    • boat principal components,

    • anchor principal components, and

    • reconstructed cube.

  • Usage Example

  • ————-

  • You can create the masks using get_annulus_segments from vip_hci.var.

  • .. code-block:: python – from vip_hci.var import get_annulus_segments

  • The function must be used as follows, where ring_out, ring_in, and

  • coro define the radius of the different annulus masks. They must have

  • the same shape as a frame of the science cube.

ones = np.ones(cube[0].shape)
boat = get_annulus_segments(ones,coro,ring_out-coro, mode="mask")[0]
anchor = get_annulus_segments(ones,ring_in,ring_out-ring_in,
                              mode="mask")[0]

Masks should be provided as ‘mask_rdi’ argument when using PCA.

res = pca(cube, angles, ref, mask_rdi=(boat, anchor), ncomp=2)

vip_hci.preproc.subsampling module

Module with pixel and frame subsampling functions.

[BRA13] (1,2)
Brandt et al. 2013
New Techniques for High-contrast Imaging with ADI: The ACORNS-ADI SEEDS Data Reduction Pipeline
The Astrophysical Journal, Volume 764, Issue 2, p. 183
vip_hci.preproc.subsampling.cube_collapse(cube, mode='median', n=50, w=None)[source]

Collapse a 3D or 4D cube into a 2D frame or 3D cube, respectively.

The mode parameter determines how the collapse should be done. It is possible to perform a trimmed mean combination of the frames, as in [BRA13]. In case of a 4D input cube, it is assumed to be an IFS dataset with the zero-th axis being the spectral dimension, and the first axis the temporal dimension.

Parameters:
  • cube (numpy ndarray) – Cube.

  • mode ({'median', 'mean', 'sum', 'max', 'trimmean', 'absmean', 'wmean'}) – Sets the way of collapsing the images in the cube. ‘wmean’ stands for weighted mean and requires weights w to be provided. ‘absmean’ stands for the mean of absolute values (potentially useful for negfc).

  • n (int, optional) – Sets the discarded values at high and low ends. When n = N is the same as taking the mean, when n = 1 is like taking the median.

  • w (1d numpy array or list, optional) – Weights to be applied for a weighted mean. Need to be provided if collapse mode is ‘wmean’.

Returns:

frame – Output array, cube combined.

Return type:

numpy ndarray

vip_hci.preproc.subsampling.cube_subsample(array, n, mode='mean', w=None, parallactic=None, verbose=True)[source]

Mean/Median combines frames in 3d or 4d cube with window n.

Parameters:
  • array (numpy ndarray) – Input 3d array, cube.

  • n (int) – Window for mean/median.

  • mode ({'mean','median', 'wmean'}, optional) – Switch for choosing mean, median or weighted average.

  • parallactic (numpy ndarray, optional) – List of corresponding parallactic angles.

  • verbose (bool optional)

Returns:

  • arr_view (numpy ndarray) – Resulting array.

  • If parallactic is provided the the new cube and angles are returned.

vip_hci.preproc.subsampling.cube_subsample_trimmean(arr, n, m)[source]

Perform a trimmed mean combination every m frames in a cube.

Details in [BRA13].

Parameters:
  • arr (numpy ndarray) – Cube.

  • n (int) – Sets the discarded values at high and low ends. When n = N is the same as taking the mean, when n = 1 is like taking the median.

  • m (int) – Window from the trimmed mean.

Returns:

arr_view – Output array, cube combined.

Return type:

numpy ndarray

Module contents

The subpackage preproc contains some useful cosmetics and pre-processing functionalities:

  • resizing frames and cubes : upscaling/pixel binning,

  • shifting frames,

  • rotating frames and cubes,

  • cropping frames and cubes,

  • removing bad pixels from frames and cubes,

  • correcting nan values from frames and cubes,

  • detecting bad frames in cubes, using:
    • pixel statistics in annulus or circular aperture,

    • ellipticity of a point like source,

    • frames correlation,

  • temporal subsampling of cubes (mean, median, trimmed mean),

  • registration (re-centering) of frames, using:
    • centroid fitting a 2d gaussian or moffat,

    • DFT upsampling or fourier cross-correlation (Guizar et al. 2008),

    • radon transform for broadband frames (Pueyo et al. 2014),

    • using satellite/waffle spots (fitting plus intersection).

  • sky subtraction (PCA method).

Astronomical calibration functionalities like flat-fielding and dark-sky subtraction, in spite of their simplicity were not included in VIP because of the heterogeneity of the datasets coming from different observatories (each having different data storage and headers). You can perform this in python in procedures of a few lines or using dedicated instrument pipelines such as esorex (ESO instruments).