vip_hci package

Subpackages

Submodules

vip_hci.hci_dataset module

Module with Dataset and Frame classes.

class vip_hci.hci_dataset.Dataset(cube, hdu=0, angles=None, wavelengths=None, fwhm=None, px_scale=None, psf=None, psfn=None, cuberef=None)[source]

Bases: vip_hci.config.utils_conf.Saveable

High-contrast imaging dataset class.

Parameters:
  • cube (str or numpy array) – 3d or 4d high-contrast image sequence. If a string is provided, cube is interpreted as the path of the FITS file containing the sequence.
  • hdu (int, optional) – If cube is a String, hdu indicates the HDU from the FITS file. By default the first HDU is used.
  • angles (list or numpy array, optional) – The vector of parallactic angles.
  • wavelengths (list or numpy array, optional) – The vector of wavelengths (to be used as scaling factors).
  • fwhm (float, optional) – The FWHM associated with this dataset (instrument dependent). Required for several methods (operations on the cube).
  • px_scale (float, optional) – The pixel scale associated with this dataset (instrument dependent).
  • psf (numpy array, optional) – The PSF template associated with this dataset.
  • psfn (numpy array, optional) – Normalized/cropped/centered version of the PSF template associated with this dataset.
  • cuberef (str or numpy array) – 3d or 4d high-contrast image sequence. To be used as a reference cube.
collapse(mode='median', n=50)[source]

Collapsing the sequence into a 2d array.

copy(deep=True, check_mem=True)[source]

Create an in-memory copy of this Dataset.

This is especially useful for keeping a backup copy of the original dataset before modifying if (e.g. with injections).

Parameters:
  • deep (bool, optional) – By default, a deep copy is created. That means every (sub)attribute is copied in memory. While this requires more memory, one can safely modify the attributes without touching the original Dataset. When deep=False, a shallow copy of the Dataset is returned instead. That means all attributes (e.g. self.cube) point back to the original object’s attributes. Pay attention when modifying such a shallow copy!
  • check_mem (bool, optional) – [deep=True] If True, verifies that the system has enough memory to store the result.
Returns:

new_dataset – (deep) copy of this Dataset.

Return type:

Dataset

crop_frames(size, xy=None, force=False)[source]

Cropping the frames of the sequence (3d or 4d cube).

Parameters:
  • size (int) – New size of the (square) frames.
  • 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.
derotate(imlib='vip-fft', interpolation='lanczos4', cxy=None, nproc=1, border_mode='constant', mask_val=nan, edge_blend=None, interp_zeros=False, ker=1)[source]

Derotating the frames of the sequence according to the parallactic angles.

Parameters:
  • imlib ({'opencv', 'skimage', 'vip-fft'}, str optional) – Library used for image transformations. Opencv is faster than ndimage or skimage.
  • interpolation (str, optional) – For ‘skimage’ 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.
  • 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 (float, optional) – 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.
drop_frames(n, m, verbose=True)[source]

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

The indices n and m are included and 1-based.

Examples

For a cube which has 5 frames numbered 1, 2, 3, 4, 5, calling ds.drop_frames(2, 4) would result in the frames 2, 3, 4 to be kept, so the first and the last frame would be discarded.

filter(method, mode, median_size=5, kernel_size=5, fwhm_size=5, btw_cutoff=0.2, btw_order=2, hann_cutoff=5, gauss_mode='conv', verbose=True)[source]

High/low pass filtering the frames of the cube.

Parameters:
  • method ({'lp', 'hp'}) – Low-pass or high-pass filtering.
  • mode (str) –

    Type of low/high-pass filtering. median

    [lp] applies a median low-pass filter to the image.
    gauss
    [lp] applies a Gaussian low-pass filter to the image.
    laplacian
    [hp] applies a Laplacian filter with kernel size defined by kernel_size using the Opencv library.
    laplacian-conv
    [hp] applies a Laplacian high-pass filter by defining a kernel (with kernel_size) and using the convolve_fft Astropy function.
    median-subt
    [hp] subtracts a median low-pass filtered version of the image.
    gauss-subt
    [hp] subtracts a Gaussian low-pass filtered version of the image.
    fourier-butter
    [hp] applies a high-pass 2D Butterworth filter in Fourier domain.
    hann
    [hp] uses a Hann window.
  • median_size (int, optional) – Size of the median box for the median or median-subt filter.
  • kernel_size (int, optional) – Size of the Laplacian kernel used in laplacian mode. It must be a positive odd integer value.
  • fwhm_size (float, optional) – Size of the Gaussian kernel used in gauss or gauss-subt mode.
  • btw_cutoff (float, optional) – Frequency cutoff for low-pass 2d Butterworth filter used in fourier-butter mode.
  • btw_order (int, optional) – Order of low-pass 2d Butterworth filter used in fourier-butter mode.
  • hann_cutoff (float) – Frequency cutoff for the hann mode.
  • gauss_mode ({'conv', 'convfft'}, str optional) – ‘conv’ uses the multidimensional gaussian filter from scipy.ndimage and ‘convfft’ uses the fft convolution with a 2d Gaussian kernel.
frame_distances(frame, region='full', dist='sad', inner_radius=None, width=None, plot=True)[source]

Calculating the frame distance/correlation with respect to a reference image.

Parameters:
  • frame (int or 2d array) – Reference frame in the cube or 2d array.
  • region ({'full', 'annulus'}, string optional) – Whether to use the full frames or a centered annulus.
  • dist ({'sad','euclidean','mse','pearson','spearman', 'ssim'}, str optional) – Which criterion to use.
  • inner_radius (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) – Whether to plot the distances or not.
frame_stats(region='circle', radius=5, xy=None, annulus_inner_radius=0, annulus_width=5, wavelength=0, plot=True)[source]

Calculating statistics on a region (circular aperture or annulus) of each image of the sequence.

Parameters:
  • region ({'circle', 'annulus'}, str optional) – Region in which basic statistics (mean, stddev, median and max) are calculated.
  • radius (int, optional) – Radius of the circular aperture.
  • xy (tuple of floats, optional) – Center of the circular aperture.
  • annulus_inner_radius (int, optional) – Inner radius of the annular region.
  • annulus_width (int, optional) – Width of the annular region.
  • wavelength (int, optional) – Index of the wavelength to be analyzed in the case of a 4d cube.
  • plot (bool, optional) – Whether to plot the frame, histograms and region.
generate_copies_with_injections(n_copies, inrad=8, outrad=12, dist_flux=('uniform', 2, 500))[source]

Create copies of this dataset, containing different random injections.

Parameters:
  • n_copies (int) – This is the number of ‘cube copies’ returned.
  • inrad,outrad (float) – Inner and outer radius of the injections. The actual injection position is chosen randomly.
  • dist_flux (tuple(‘method’, *params)) –

    Tuple describing the flux selection. Method can be a function, the *params are passed to it. Method can also be a string, for a pre-defined random function:

    ('skewnormal', skew, mean, var)
    uses scipy.stats.skewnorm.rvs
    ('uniform', low, high)
    uses np.random.uniform
    ('normal', loc, scale)
    uses np.random.normal
Yields:

fake_dataset (Dataset) – Copy of the original Dataset, with injected companions.

get_nbytes()[source]

Return the total number of bytes the Dataset consumes.

inject_companions(flux, rad_dists, n_branches=1, theta=0, imlib='vip-fft', interpolation='lanczos4', full_output=False, verbose=True)[source]

Injection of fake companions in 3d or 4d cubes.

Parameters:
  • flux (float or list) – Factor for controlling the brightness of the fake companions.
  • rad_dists (float, list or array 1d) – Vector of radial distances of fake companions in pixels.
  • n_branches (int, optional) – Number of azimutal branches.
  • theta (float, optional) – Angle in degrees for rotating the position of the first branch that by default is located at zero degrees. Theta counts counterclockwise from the positive x axis.
  • imlib ({'opencv', 'ndimage-fourier', 'ndimage-interp', 'vip-fft'}, str opt) – Library or method used for performing the image shift. ‘ndimage-fourier’, does a fourier shift operation and preserves better the pixel values (therefore the flux and photometry). Interpolation based shift (‘opencv’ and ‘ndimage-interp’) is faster than the fourier shift. ‘opencv’ is recommended when speed is critical.
  • interpolation ({'bicubic', 'bilinear', 'nearneig'}, optional) – Only used in case of imlib is set to ‘opencv’ or ‘ndimage-interp’, where the images are shifted via interpolation. For ‘ndimage-interp’ 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.
  • full_output (bool, optional) – Return the coordinates of the injected companions.
  • verbose (bool, optional) – If True prints out additional information.
Returns:

yx – [full_output=True] Pixel coordinates of the injections in the first frame (and first wavelength for 4D cubes). These are only the new injections - all injections (from multiple calls to this function) are stored in self.injections_yx.

Return type:

list of tuple(y,x)

load_angles(angles, hdu=0)[source]

Loads the PA vector from a FITS file. It is possible to specify the HDU.

Parameters:
  • angles (str or 1d numpy ndarray) – List or vector with the parallactic angles.
  • hdu (int, optional) – If angles is a String, hdu indicates the HDU from the FITS file. By default the first HDU is used.
load_wavelengths(wavelengths, hdu=0)[source]

Loads the scaling factors vector from a FITS file. It is possible to specify the HDU.

Parameters:
  • wavelengths (str or 1d numpy ndarray) – List or vector with the wavelengths.
  • hdu (int, optional) – If wavelengths is a String, hdu indicates the HDU from the FITS file. By default the first HDU is used.
mask_center(radius, fillwith=0, mode='in')[source]

Masking the values inside/outside a centered circular aperture.

Parameters:
  • radius (int) – Radius of the circular aperture.
  • fillwith (int, float or np.nan, optional) – Value to put instead of the masked out pixels.
  • mode ({'in', 'out'}, optional) – When set to ‘in’ then the pixels inside the radius are set to fillwith. When set to ‘out’ the pixels outside the circular mask are set to fillwith.
normalize_psf(fit_fwhm=True, size=None, threshold=None, mask_core=None, model='gauss', imlib='vip-fft', interpolation='lanczos4', force_odd=True, verbose=True)[source]

Normalizes a PSF (2d or 3d array), to have the flux in a 1xFWHM aperture equal to one. It also allows to crop the array and center the PSF at the center of the frame(s).

Parameters:
  • fit_fwhm (bool, optional) – Whether to fit a model to estimate the FWHM instead of using the self.fwhm attribute.
  • size (int or None, optional) – If int it will correspond to the size of the squared subimage to be cropped form the psf array.
  • threshold (None of float, optional) – Sets to zero small values, trying to leave only the core of the PSF.
  • mask_core (None of float, optional) – Sets the radius of a circular aperture for the core of the PSF, everything else will be set to zero.
  • imlib ({'opencv', 'ndimage-fourier', 'ndimage-interp', 'vip-fft'}, str opt) – Library or method used for performing the image shift. ‘ndimage-fourier’, does a fourier shift operation and preserves better the pixel values (therefore the flux and photometry). Interpolation based shift (‘opencv’ and ‘ndimage-interp’) is faster than the fourier shift. ‘opencv’ is recommended when speed is critical.
  • interpolation ({'bicubic', 'bilinear', 'nearneig'}, optional) – Only used in case of imlib is set to ‘opencv’ or ‘ndimage-interp’, where the images are shifted via interpolation. For ‘ndimage-interp’ 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.
  • force_odd (str, optional) – If True the resulting array will have odd size (and the PSF will be placed at its center). If False, and the frame size is even, then the PSF will be put at the center of an even-sized frame.
  • verbose (bool, optional) – If True intermediate results are printed out.
plot(**kwargs)[source]

Plotting the frames of a 3D or 4d cube.

Parameters:**kwargs (dict, optional) – Parameters passed to the function plot_cubes of the package HCIplot.
recenter(method='2dfit', xy=None, subi_size=5, model='gauss', nproc=1, imlib='vip-fft', interpolation='lanczos4', offset=None, negative=False, threshold=False, save_shifts=False, cy_1=None, cx_1=None, upsample_factor=100, alignment_iter=5, gamma=1, min_spat_freq=0.5, max_spat_freq=3, recenter_median=False, sigfactor=6, cropsize=101, hsize=0.4, step=0.01, mask_center=None, verbose=True, debug=False, plot=True)[source]

Frame to frame recentering.

Parameters:
  • method ({'2dfit', 'dftups', 'dftupspeckles', 'satspots', 'radon'}, optional) – Recentering method.
  • xy (tuple or ints or tuple of 4 tuples of ints, optional) – For the 2dfitting, xy are the coordinates of the center of the subimage (wrt the original frame). For the satellite spots method, it is a 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 the square subimage sides in pixels.
  • model (str, optional) – [method=2dfit] Sets the type of fit to be used. ‘gauss’ for a 2d Gaussian fit and ‘moff’ for a 2d Moffat 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 ({'opencv', 'ndimage-fourier', 'ndimage-interp', 'vip-fft'}, str opt) – Library or method used for performing the image shift. ‘ndimage-fourier’, does a fourier shift operation and preserves better the pixel values (therefore the flux and photometry). Interpolation based shift (‘opencv’ and ‘ndimage-interp’) is faster than the fourier shift. ‘opencv’ is recommended when speed is critical.
  • interpolation ({'bicubic', 'bilinear', 'nearneig'}, optional) – Only used in case of imlib is set to ‘opencv’ or ‘ndimage-interp’, where the images are shifted via interpolation. For ‘ndimage-interp’ 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.
  • offset (tuple of floats, optional) – [method=2dfit] 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) – [method=2dfit/dftups/dftupspeckles] If True a negative 2d Gaussian/Moffat fit is performed.
  • threshold (bool, optional) – [method=2dfit] If True the background pixels (estimated using sigma clipped statistics) will be replaced by small random Gaussian noise.
  • save_shifts (bool, optional) – [method=2dfit/dftups] Whether to save the shifts to a file in disk.
  • cx_1 (cy_1,) – [method=dftups] Coordinates of the center of the subimage for fitting a 2d Gaussian and centroiding the 1st frame.
  • upsample_factor (int, optional) – [method=dftups] Upsampling factor (default 100). Images will be registered to within 1/upsample_factor of a pixel.
  • alignment_iter (int, optional) – [method=dftupspeckles] Number of alignment iterations (recomputes median after each iteration).
  • gamma (int, optional) – [method=dftupspeckles] Applies a gamma correction to emphasize speckles (useful for faint stars).
  • min_spat_freq (float, optional) – [method=dftupspeckles] Spatial frequency for high pass filter.
  • max_spat_freq (float, optional) – [method=dftupspeckles] Spatial frequency for low pass filter.
  • recenter_median (bool, optional) – [method=dftupspeckles] Recenter the frames at each iteration based on the gaussian fit.
  • sigfactor (int, optional) – [method=satspots] 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.
  • cropsize (odd int, optional) – [method=radon] Size in pixels of the cropped central area of the input array that will be used. It should be large enough to contain the satellite spots.
  • hsize (float, optional) – [method=radon] Size of the box for the grid search. The frame is shifted to each direction from the center in a hsize length with a given step.
  • step (float, optional) – [method=radon] The step of the coordinates change.
  • mask_center (None or int, optional) – [method=radon] 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.
  • verbose (bool, optional) – Whether to print to stdout the timing and aditional info.
  • debug (bool, optional) – If True debug information is printed and plotted.
  • plot (bool, optional) – Whether to plot the shifts.
remove_badframes(method='corr', frame_ref=None, crop_size=30, dist='pearson', percentile=20, stat_region='annulus', inner_radius=10, width=10, top_sigma=1.0, low_sigma=1.0, window=None, roundlo=-0.2, roundhi=0.2, lambda_ref=0, plot=True, verbose=True)[source]

Find outlying/bad frames and slice the cube accordingly.

Besides modifying self.cube and self.angles, also sets a self.good_indices which contain the indices of the angles which were kept.

Parameters:
  • method ({'corr', 'pxstats', 'ellip'}, optional) – Method which is used to determine bad frames. Refer to the preproc.badframes submodule for explanation of the different methods.
  • frame_ref (int, 2d array or None, optional) – [method=corr] Index of the frame that will be used as a reference or 2d reference array.
  • crop_size (int, optional) – [method=corr] Size in pixels of the square subframe to be analyzed.
  • dist ({'sad','euclidean','mse','pearson','spearman'}, optional) – [method=corr] One of the similarity or dissimilarity measures from function vip_hci.stats.distances.cube_distance().
  • percentile (float, optional) – [method=corr] The percentage of frames that will be discarded [0..100].
  • stat_region ({'annulus', 'circle'}, optional) – [method=pxstats] Whether to take the statistics from a circle or an annulus.
  • inner_radius (int, optional) – [method=pxstats] If stat_region is ‘annulus’ then ‘in_radius’ is the inner radius of the annular region. If stat_region is ‘circle’ then ‘in_radius’ is the radius of the aperture.
  • width (int, optional) – [method=pxstats] Size of the annulus. Ignored if mode is ‘circle’.
  • top_sigma (int, optional) – [method=pxstats] Top boundary for rejection.
  • low_sigma (int, optional) – [method=pxstats] Lower boundary for rejection.
  • window (int, optional) – [method=pxstats] Window for smoothing the median and getting the rejection statistic.
  • roundlo,roundhi (float, optional) – [method=ellip] : Lower and higher bounds for the ellipticipy.
  • lambda_ref (int, optional) – [4D cube] Which wavelength to consider when determining bad frames on a 4D cube.
  • plot (bool, optional) – If true it plots the mean fluctuation as a function of the frames and the boundaries.
  • verbose (bool, optional) – Show debug output.
rescale(scale, imlib='ndimage', interpolation='bicubic', verbose=True)[source]

Resampling the pixels (upscaling or downscaling the frames).

Parameters:
  • 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 ({'ndimage', 'opencv', 'vip-fft'}, str optional) – Library used for image transformations. ndimage is the default.
  • 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 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.
  • verbose (bool, optional) – Whether to print out additional info such as the new cube shape.
subsample(window, mode='mean')[source]

Temporally sub-sampling the sequence (3d or 4d cube).

Parameters:
  • window (int) – Window for mean/median.
  • mode ({'mean', 'median'}, optional) – Switch for choosing mean or median.
class vip_hci.hci_dataset.Frame(data, hdu=0, fwhm=None)[source]

Bases: object

High-contrast imaging frame (2d array).

Parameters:
  • data (numpy ndarray) – 2d array.
  • hdu (int, optional) – If cube is a String, hdu indicates the HDU from the FITS file. By default the first HDU is used.
  • fwhm (float, optional) – The FWHM associated with this dataset (instrument dependent). Required for several methods (operations on the cube).
crop(size, xy=None, force=False)[source]

Cropping the frame.

Parameters:
  • 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 True this condition can be avoided.
detect_blobs(psf, bkg_sigma=1, method='lpeaks', matched_filter=False, mask=True, snr_thresh=5, plot=True, debug=False, verbose=False, save_plot=None, plot_title=None, angscale=False)[source]

Detecting blobs on the 2d array.

filter(method, mode, median_size=5, kernel_size=5, fwhm_size=5, btw_cutoff=0.2, btw_order=2, hann_cutoff=5, gauss_mode='conv')[source]

High/low pass filtering the frames of the image.

Parameters:
  • method ({'lp', 'hp'}) – Low-pass or high-pass filtering.
  • mode (str) –

    Type of low/high-pass filtering. median

    [lp] applies a median low-pass filter to the image.
    gauss
    [lp] applies a Gaussian low-pass filter to the image.
    laplacian
    [hp] applies a Laplacian filter with kernel size defined by kernel_size using the Opencv library.
    laplacian-conv
    [hp] applies a Laplacian high-pass filter by defining a kernel (with kernel_size) and using the convolve_fft Astropy function.
    median-subt
    [hp] subtracts a median low-pass filtered version of the image.
    gauss-subt
    [hp] subtracts a Gaussian low-pass filtered version of the image.
    fourier-butter
    [hp] applies a high-pass 2D Butterworth filter in Fourier domain.
    hann
    [hp] uses a Hann window.
  • median_size (int, optional) – Size of the median box for the median or median-subt filter.
  • kernel_size (int, optional) – Size of the Laplacian kernel used in laplacian mode. It must be a positive odd integer value.
  • fwhm_size (float, optional) – Size of the Gaussian kernel used in gauss or gauss-subt mode.
  • btw_cutoff (float, optional) – Frequency cutoff for low-pass 2d Butterworth filter used in fourier-butter mode.
  • btw_order (int, optional) – Order of low-pass 2d Butterworth filter used in fourier-butter mode.
  • hann_cutoff (float) – Frequency cutoff for the hann mode.
  • gauss_mode ({'conv', 'convfft'}, str optional) – ‘conv’ uses the multidimensional gaussian filter from scipy.ndimage and ‘convfft’ uses the fft convolution with a 2d Gaussian kernel.
get_center(verbose=True)[source]

Getting the coordinates of the center of the image.

Parameters:verbose (bool optional) – If True the center coordinates are printed out.
plot(**kwargs)[source]

Plotting the 2d array.

Parameters:**kwargs (dict, optional) – Parameters passed to the function plot_frames of the package HCIplot.
radial_profile(sep=1)[source]

Calculates the average radial profile of an image.

Parameters:sep (int, optional) – The average radial profile is recorded every sep pixels.
recenter(method='satspots', xy=None, subi_size=19, sigfactor=6, imlib='vip-fft', interpolation='lanczos4', debug=False, verbose=True)[source]

Recentering the frame using the satellite spots or a radon transform.

Parameters:
  • method ({'satspots', 'radon'}, str optional) – Method for recentering the frame.
  • xy (tuple, optional) – Tuple with coordinates X,Y of the 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.
rescale(scale, imlib='vip-fft', interpolation='bicubic', verbose=True)[source]

Resampling the image (upscaling or downscaling).

Parameters:
  • 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 ({'ndimage', 'opencv', 'vip-fft'}, str optional) – Library used for image transformations. ndimage is the default.
  • 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 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.
  • verbose (bool, optional) – Whether to print out additional info such as the new image shape.
rotate(angle, imlib='vip-fft', interpolation='lanczos4', cxy=None)[source]

Rotating the image by a given angle.

Parameters:
  • imlib ({'opencv', 'skimage', 'vip-fft'}, str optional) – Library used for image transformations. Opencv is faster than ndimage or skimage.
  • interpolation (str, optional) – For ‘skimage’ 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.
  • 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.
shift(shift_y, shift_x, imlib='vip-fft', interpolation='lanczos4')[source]

Shifting the image.

Parameters:
  • shift_x (shift_y,) – Shifts in x and y directions.
  • imlib ({'opencv', 'ndimage-fourier', 'ndimage-interp', 'vip-fft'}, str opt) – Library or method used for performing the image shift. ‘ndimage-fourier’, does a fourier shift operation and preserves better the pixel values (therefore the flux and photometry). Interpolation based shift (‘opencv’ and ‘ndimage-interp’) is faster than the fourier shift. ‘opencv’ is recommended when speed is critical.
  • interpolation ({'bicubic', 'bilinear', 'nearneig'}, optional) – Only used in case of imlib is set to ‘opencv’ or ‘ndimage-interp’, where the images are shifted via interpolation. For ‘ndimage-interp’ 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.
snr(source_xy, plot=False, verbose=True)[source]

Calculating the S/N for a test resolution element source_xy.

Parameters:
  • source_xy (tuple of floats) – X and Y coordinates of the planet or test speckle.
  • plot (bool, optional) – Plots the frame and the apertures considered for clarity.
  • verbose (bool, optional) – Chooses whether to print some output or not.
Returns:

snr_val – Value of the S/N for source_xy.

Return type:

float

stats(region='circle', radius=5, xy=None, annulus_inner_radius=0, annulus_width=5, source_xy=None, verbose=True, plot=True)[source]

Calculating statistics on the image, both in the full-frame and in a region (circular aperture or annulus). Also, the S/N of the either source_xy or the max pixel is calculated.

Parameters:
  • region ({'circle', 'annulus'}, str optional) – Region in which basic statistics (mean, stddev, median and max) are calculated.
  • radius (int, optional) – Radius of the circular aperture.
  • xy (tuple of floats, optional) – Center of the circular aperture.
  • annulus_inner_radius (int, optional) – Inner radius of the annular region.
  • annulus_width (int, optional) – Width of the annular region.
  • source_xy (tuple of floats, optional) – Coordinates for which the S/N information will be obtained. If None, the S/N is estimated for the pixel with the maximum value.
  • verbose (bool, optional) – Whether to print out the values of the calculated statistics.
  • plot (bool, optional) – Whether to plot the frame, histograms and region.

vip_hci.hci_postproc module

Module with the HCI<post-processing algorithms> classes.

class vip_hci.hci_postproc.HCIMedianSub(dataset=None, mode='fullfr', radius_int=0, asize=1, delta_rot=1, delta_sep=(0.2, 1), nframes=4, imlib='vip-fft', interpolation='lanczos4', collapse='median', verbose=True)[source]

Bases: vip_hci.hci_postproc.HCIPostProcAlgo

HCI median subtraction algorithm.

Parameters:
  • mode ({"fullfr","annular"}, str optional) – In “simple” mode only the median frame is subtracted, in “annular” mode also the 4 closest frames given a PA threshold (annulus-wise) are subtracted.
  • radius_int (int, optional) – The radius of the innermost annulus. By default is 0, if >0 then the central circular area is discarded.
  • asize (int, optional) – The size of the annuli, in FWHM. Default is 2.
  • delta_rot (int, optional) – Factor for increasing the parallactic angle threshold, expressed in FWHM. Default is 1 (excludes 1 FHWM on each side of the considered frame).
  • delta_sep (float or tuple of floats, optional) – The threshold separation in terms of the mean FWHM (for ADI+mSDI data). If a tuple of two values is provided, they are used as the lower and upper intervals for the threshold (grows as a function of the separation).
  • nframes (even int, optional) – Number of frames to be used for building the optimized reference PSF when working in annular mode.
  • imlib ({'opencv', 'skimage', 'vip-fft'}, str optional) – Library used for image transformations. Opencv is faster than ndimage or skimage.
  • interpolation (str, optional) – For ‘skimage’ 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.
  • verbose (bool, optional) – Show more output.
run(dataset=None, nproc=1, verbose=True)[source]

Running the HCI median subtraction algorithm for model PSF subtraction.

Parameters:
  • dataset (HCIDataset object) – An HCIDataset object to be processed.
  • full_output (bool, optional) – Whether to return the final median combined image only or with other intermediate arrays.
  • nproc (int, optional) – Number of processes for parallel computing. Defaults to single-core processing.
  • verbose (bool, optional) – If True prints to stdout intermediate info.
class vip_hci.hci_postproc.HCIPca(dataset=None, ncomp=1, ncomp2=1, svd_mode='lapack', scaling=None, adimsdi='double', mask_center_px=None, source_xy=None, delta_rot=1, imlib='vip-fft', interpolation='lanczos4', collapse='median', check_mem=True, crop_ifs=True, verbose=True)[source]

Bases: vip_hci.hci_postproc.HCIPostProcAlgo

HCI PCA algorithm.

Parameters:
  • dataset (HCIDataset object, optional) – An HCIDataset object to be processed. Can also be passed to .run().
  • ncomp (int, optional) – How many PCs are used as a lower-dimensional subspace to project the target frames. For an ADI cube, ncomp is the number of PCs extracted from cube. For the RDI case, when cube and cube_ref are provided, ncomp is the number of PCs obtained from cube_ref. For an ADI+mSDI cube (e.g. SPHERE/IFS), if adimsdi is double then ncomp is the number of PCs obtained from each multi-spectral frame (if ncomp is None then this stage will be skipped and the spectral channels will be combined without subtraction). If adimsdi is single, then ncomp is the number of PCs obtained from the whole set of frames (n_channels * n_adiframes).
  • ncomp2 (int, optional) – Only used for ADI+mSDI cubes, when adimsdi is set to double. ncomp2 sets the number of PCs used in the second PCA stage (ADI fashion, using the residuals of the first stage). If None then the second PCA stage is skipped and the residuals are de-rotated and combined.
  • svd_mode ({'lapack', 'arpack', 'eigen', 'randsvd', 'cupy', 'eigencupy', 'randcupy'}, str) – Switch for the SVD method/library to be used. lapack uses the LAPACK linear algebra library through Numpy and it is the most conventional way of computing the SVD (deterministic result computed on CPU). arpack uses the ARPACK Fortran libraries accessible through Scipy (computation on CPU). eigen computes the singular vectors through the eigendecomposition of the covariance M.M’ (computation on CPU). randsvd uses the randomized_svd algorithm implemented in Sklearn (computation on CPU). cupy uses the Cupy library for GPU computation of the SVD as in the LAPACK version. eigencupy offers the same method as with the eigen option but on GPU (through Cupy). randcupy is an adaptation of the randomized_svd algorithm, where all the computations are done on a GPU.
  • scaling ({None, 'temp-mean', 'spat-mean', 'temp-standard', 'spat-standard'}) – With None, no scaling is performed on the input data before SVD. With “temp-mean” then temporal px-wise mean subtraction is done, with “spat-mean” then the spatial mean is subtracted, with “temp-standard” temporal mean centering plus scaling to unit variance is done and with “spat-standard” spatial mean centering plus scaling to unit variance is performed.
  • adimsdi ({'double', 'single'}, str optional) – In the case cube is a 4d array, adimsdi determines whether a single or double pass PCA is going to be computed. In the single case, the multi-spectral frames are rescaled wrt the largest wavelength to align the speckles and all the frames are processed with a single PCA low-rank approximation. In the double case, a firt stage is run on the rescaled spectral frames, and a second PCA frame is run on the residuals in an ADI fashion.
  • mask_center_px (None or int) – If None, no masking is done. If an integer > 1 then this value is the radius of the circular mask.
  • source_xy (tuple of int, optional) – For ADI PCA, this triggers a frame rejection in the PCA library. source_xy are the coordinates X,Y of the center of the annulus where the PA criterion will be used to reject frames from the library.
  • delta_rot (int, optional) – Factor for tunning the parallactic angle threshold, expressed in FWHM. Default is 1 (excludes 1xFHWM on each side of the considered frame).
  • imlib ({'opencv', 'skimage', 'vip-fft'}, str optional) – Library used for image transformations. Opencv is faster than ndimage or skimage.
  • interpolation (str, optional) – For ‘skimage’ 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.
  • check_mem (bool, optional) – If True, it check that the input cube(s) are smaller than the available system memory.
run(dataset=None, nproc=1, verbose=True, debug=False)[source]

Run the HCI PCA algorithm for model PSF subtraction.

Note

creates/sets the self.frame_final attribute, and depending on the input parameters:

3D case:
cube_reconstructed cube_residuals cube_residuals_der
3D case, source_xy is None:
cube_residuals pcs
4D case, adimsdi=”double”:
cube_residuals_per_channel cube_residuals_per_channel_der
4D case, adimsdi=”single”:
cube_residuals cube_residuals_resc
Parameters:
  • dataset (HCIDataset, optional) – Dataset to process. If not provided, self.dataset is used (as set when initializing this object).
  • nproc (int, optional) – (not used) Note that HCIPca always works in single-processing mode.
  • verbose (bool, optional) – Show more output.
  • debug (bool, optional) – Whether to print debug information or not.
class vip_hci.hci_postproc.HCILoci(dataset=None, scale_list=None, metric='manhattan', dist_threshold=90, delta_rot=0.5, delta_sep=(0.1, 1), radius_int=0, asize=4, n_segments=4, solver='lstsq', tol=0.001, optim_scale_fact=1, adimsdi='skipadi', imlib='vip-fft', interpolation='lanczos4', collapse='median', verbose=True)[source]

Bases: vip_hci.hci_postproc.HCIPostProcAlgo

HCI LOCI algorithm.

run(dataset=None, nproc=1, verbose=True)[source]

Run the HCI LOCI algorithm for model PSF subtraction.

class vip_hci.hci_postproc.HCILLSG(dataset=None, rank=10, thresh=1, max_iter=10, low_rank_ref=False, low_rank_mode='svd', auto_rank_mode='noise', residuals_tol=0.1, cevr=0.9, thresh_mode='soft', nproc=1, asize=None, n_segments=4, azimuth_overlap=None, radius_int=None, random_seed=None, imlib='vip-fft', interpolation='lanczos4', high_pass=None, collapse='median', verbose=True)[source]

Bases: vip_hci.hci_postproc.HCIPostProcAlgo

HCI LLSG algorithm.

run(dataset=None, nproc=1, verbose=True)[source]

Run the HCI LLSG algorithm for model PSF subtraction.

class vip_hci.hci_postproc.HCIAndromeda(dataset=None, oversampling_fact=0.5, filtering_fraction=0.25, min_sep=0.5, annuli_width=1.0, roa=2.0, opt_method='lsq', nsmooth_snr=18, iwa=None, owa=None, precision=50, fast=False, homogeneous_variance=True, ditimg=1.0, ditpsf=None, tnd=1.0, total=False, multiply_gamma=True, verbose=True)[source]

Bases: vip_hci.hci_postproc.HCIPostProcAlgo

HCI ANDROMEDA algorithm.

Parameters:
  • dataset (HCIDataset object, optional) – An HCIDataset object to be processed. Can also be passed to .run().
  • oversampling_fact (float, optional) – Oversampling factor for the wavelength corresponding to the filter used for obtaining cube (defined as the ratio between the wavelength of the filter and the Shannon wavelength).
  • filtering_fraction (float, optional) – Strength of the high-pass filter. If set to 1, no high-pass filter is used.
  • min_sep (float, optional) – Angular separation is assured to be above min_sep*lambda/D.
  • annuli_width (float, optional) – Annuli width on which the subtraction are performed. The same for all annuli.
  • roa (float, optional) – Ratio of the optimization area. The optimization annulus area is defined by roa * annuli_width.
  • opt_method ({'no', 'total', 'lsq', 'robust'}, optional) – Method used to balance for the flux difference that exists between the two subtracted annuli in an optimal way during ADI.
  • nsmooth_snr (int, optional) – Number of pixels over which the radial robust standard deviation profile of the SNR map is smoothed to provide a global trend for the SNR map normalization. For nsmooth_snr=0 the SNR map normalization is disabled, and the positivity constraint is applied when calculating the flux.
  • iwa (float, optional) – Inner working angle / inner radius of the first annulus taken into account, expressed in $lambda/D$.
  • precision (int, optional) – Number of shifts applied to the PSF. Passed to calc_psf_shift_subpix , which then creates a 4D cube with shape (precision+1, precision+1, N, N).
  • homogeneous_variance (bool, optional) – If set, variance is treated as homogeneous and is calculated as a mean of variance in each position through time.
  • multiply_gamma (bool, optional) – Use gamma for signature computation too.
  • verbose (bool, optional) – Print some parameter values for control.
make_snr_map(*args, **kwargs)[source]

Does nothing. For Andromeda, snr_map is calculated by run().

Note

The @calculates decorator is not present in this function definition, so self._get_calculations does not mark snr_map as created by this function.

run(dataset=None, nproc=1, verbose=True)[source]

Run the ANDROMEDA algorithm for model PSF subtraction.

Parameters:
  • dataset (HCIDataset, optional) – Dataset to process. If not provided, self.dataset is used (as set when initializing this object).
  • nproc (int, optional) – Number of processes to use.
  • verbose (bool, optional) – Print some parameter values for control.

vip_hci.vip_ds9 module

Module with a class for creating a DS9 window through pyds9.

Module contents