vip_hci.psfsub package

Submodules

vip_hci.psfsub.framediff module

Module with a frame differencing algorithm for ADI post-processing.

class vip_hci.psfsub.framediff.FRAME_DIFF_Params(cube: ndarray | None = None, angle_list: ndarray | None = None, fwhm: float = 4, metric: Enum = Metric.MANHATTAN, dist_threshold: int = 50, n_similar: int | None = None, delta_rot: int = 0.5, radius_int: int = 2, asize: int = 4, ncomp: int | None = None, imlib: Enum = Imlib.VIPFFT, interpolation: Enum = Interpolation.LANCZOS4, collapse: Enum = Collapse.MEDIAN, nproc: int = 1, verbose: bool = True, debug: bool = False, full_output: bool = False)[source]

Bases: object

Set of parameters for the frame differencing module.

See the function frame_diff below for the documentation.

angle_list: ndarray = None
asize: int = 4
collapse: Enum = 'median'
cube: ndarray = None
debug: bool = False
delta_rot: int = 0.5
dist_threshold: int = 50
full_output: bool = False
fwhm: float = 4
imlib: Enum = 'vip-fft'
interpolation: Enum = 'lanczos4'
metric: Enum = 'manhattan'
n_similar: int = None
ncomp: int = None
nproc: int = 1
radius_int: int = 2
verbose: bool = True
vip_hci.psfsub.framediff.frame_diff(*all_args: List, **all_kwargs: dict)[source]

Run the frame differencing algorithm.

It uses vector distance (depending on metric), using separately the pixels from different annuli of asize width, to create pairs of most similar images. Then it performs pair-wise subtraction and combines the residuals.

Parameters:
  • all_args (list, optional) – Positionnal arguments for the frame diff algorithm. Full list of parameters below.

  • all_kwargs (dictionary, optional) – Mix of keyword arguments that can initialize a FrameDiffParams and the optional ‘rot_options’ dictionnary, with keyword values for “border_mode”, “mask_val”, “edge_blend”, “interp_zeros”, “ker” (see documentation of vip_hci.preproc.frame_rotate). Can also contain a FrameDiffParams named as algo_params.

  • parameters (Frame differencing)

  • ----------

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

  • angle_list (numpy ndarray, 1d) – Corresponding parallactic angle for each frame.

  • fwhm (float, optional) – Known size of the FWHM in pixels to be used. Default is 4.

  • metric (str, optional) – Distance metric to be used (‘cityblock’, ‘cosine’, ‘euclidean’, ‘l1’, ‘l2’, ‘manhattan’, ‘correlation’, etc). It uses the scikit-learn function sklearn.metrics.pairwise.pairwise_distances (check its documentation).

  • dist_threshold (int) – Indices with a distance larger than dist_threshold percentile will initially discarded.

  • n_similar (None or int, optional) – If a postive integer value is given, then a median combination of n_similar frames will be used instead of the most similar one.

  • delta_rot (int) – Minimum parallactic angle distance between the pairs.

  • 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 pixels.

  • ncomp (None or int, optional) – If a positive integer value is given, then the annulus-wise PCA low-rank approximation with ncomp principal components will be subtracted. The pairwise subtraction will be performed on these residuals.

  • nproc (None or int, optional) – Number of processes for parallel computing. If None the number of processes will be set to cpu_count()/2. By default the algorithm works in single-process mode.

  • imlib (str, opt) – See description in vip.preproc.frame_rotate()

  • interpolation (str, opt) – See description in vip.preproc.frame_rotate()

  • collapse (str, opt) – What to do with derotated residual cube? See options of vip.preproc.cube_collapse()

  • verbose (bool, optional) – If True prints info to stdout.

  • debug (bool, optional) – If True the distance matrices will be plotted and additional information will be given.

Returns:

final_frame – Median combination of the de-rotated cube.

Return type:

numpy ndarray, 2d

vip_hci.psfsub.llsg module

Module containing the Local Low-rank plus Sparse plus Gaussian-noise decomposition algorithm for ADI data.

[GOM16]
Gomez-Gonzalez et al. 2016
Low-rank plus sparse decomposition for exoplanet detection in direct-imaging ADI sequences. The LLSG algorithm
Astronomy & Astrophysics, Volume 589, Issue 1, p. 54
class vip_hci.psfsub.llsg.LLSG_Params(cube: ndarray | None = None, angle_list: ndarray | None = None, fwhm: float | None = None, rank: int = 10, thresh: float = 1, max_iter: int = 10, low_rank_ref: bool = False, low_rank_mode: Enum = LowRankMode.SVD, auto_rank_mode: Enum = AutoRankMode.NOISE, residuals_tol: float = 0.1, cevr: float = 0.9, thresh_mode: Enum = ThreshMode.SOFT, nproc: int = 1, asize: int | None = None, n_segments: int = 4, azimuth_overlap: int | None = None, radius_int: int | None = None, random_seed: int | None = None, high_pass: int | None = None, collapse: Enum = Collapse.MEDIAN, full_output: bool = False, verbose: bool = True, debug: bool = False)[source]

Bases: object

Set of parameters for the LLSG algorithm.

See function llsg below for the documentation.

angle_list: ndarray = None
asize: int = None
auto_rank_mode: Enum = 'noise'
azimuth_overlap: int = None
cevr: float = 0.9
collapse: Enum = 'median'
cube: ndarray = None
debug: bool = False
full_output: bool = False
fwhm: float = None
high_pass: int = None
low_rank_mode: Enum = 'svd'
low_rank_ref: bool = False
max_iter: int = 10
n_segments: int = 4
nproc: int = 1
radius_int: int = None
random_seed: int = None
rank: int = 10
residuals_tol: float = 0.1
thresh: float = 1
thresh_mode: Enum = 'soft'
verbose: bool = True
vip_hci.psfsub.llsg.llsg(*all_args: List, **all_kwargs: dict)[source]

Local Low-rank plus Sparse plus Gaussian-noise decomposition (LLSG) as described in [GOM16]. This first version of our algorithm aims at decomposing ADI cubes into three terms L+S+G (low-rank, sparse and Gaussian noise). Separating the noise from the S component (where the moving planet should stay) allow us to increase the SNR of potential planets.

The three tunable parameters are the rank or expected rank of the L component, the thresh or threshold for encouraging sparsity in the S component and max_iter which sets the number of iterations. The rest of parameters can be tuned at the users own risk (do it if you know what you’re doing).

Parameters:
  • all_args (list, optional) – Positionnal arguments for the LLSG algorithm. Full list of parameters below.

  • all_kwargs (dictionary, optional) – Mix of keyword arguments that can initialize a LLSGParams and the optional ‘rot_options’ dictionnary, with keyword values for “border_mode”, “mask_val”, “edge_blend”, “interp_zeros”, “ker” (see documentation of vip_hci.preproc.frame_rotate). Can also contain a LLSGParams named as algo_params.

  • parameters (LLSG)

  • ----------

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

  • angle_list (numpy ndarray, 1d) – Corresponding parallactic angle for each frame.

  • fwhm (float) – Known size of the FWHM in pixels to be used.

  • rank (int, optional) – Expected rank of the L component.

  • thresh (float, optional) – Factor that scales the thresholding step in the algorithm.

  • max_iter (int, optional) – Sets the number of iterations.

  • low_rank_ref – If True the first estimation of the L component is obtained from the remaining segments in the same annulus.

  • low_rank_mode (Enum, see vip_hci.config.paramenum.LowRankMode) – Sets the method of solving the L update.

  • auto_rank_mode (Enum, see vip_hci.config.paramenum.AutoRankMode) – If rank is None, then auto_rank_mode sets the way that the rank is determined: the noise minimization or the cumulative explained variance ratio (when ‘svd’ is used).

  • residuals_tol (float, optional) – The value of the noise decay to be used when rank is None and auto_rank_mode is set to noise.

  • cevr (float, optional) – Float value in the range [0,1] for selecting the cumulative explained variance ratio to choose the rank automatically (if rank is None).

  • thresh_mode (Enum, see vip_hci.config.paramenum.ThreshMode) – Sets the type of thresholding.

  • nproc (None or int, optional) – Number of processes for parallel computing. If None the number of processes will be set to cpu_count()/2. By default the algorithm works in single-process mode.

  • asize (int or None, optional) – If asize is None then each annulus will have a width of 2*asize. If an integer then it is the width in pixels of each annulus.

  • n_segments (int or list of ints, optional) – The number of segments for each annulus. When a single integer is given it is used for all annuli.

  • azimuth_overlap (int or None, optional) – Sets the amount of azimuthal averaging.

  • radius_int (int, optional) – The radius of the innermost annulus. By default is 0, if >0 then the central circular area is discarded.

  • random_seed (int or None, optional) – Controls the seed for the Pseudo Random Number generator.

  • high_pass (odd int or None, optional) – If set to an odd integer <=7, a high-pass filter is applied to the frames. The vip_hci.var.frame_filter_highpass is applied twice, first with the mode median-subt and a large window, and then with laplacian-conv and a kernel size equal to high_pass. 5 is an optimal value when fwhm is ~4.

  • collapse (Enum, see vip_hci.config.paramenum.Collapse) – Sets the way of collapsing the frames for producing a final image.

  • full_output (bool, optional) – Whether to return the final median combined image only or with other intermediate arrays.

  • verbose (bool, optional) – If True prints to stdout intermediate info.

  • debug (bool, optional) – Whether to output some intermediate information.

Returns:

  • frame_s (numpy ndarray, 2d) – Final frame (from the S component) after rotation and median-combination.

  • If full_output is True, the following intermediate arrays are returned

  • list_l_array_der, list_s_array_der, list_g_array_der, frame_l, frame_s,

  • frame_g

vip_hci.psfsub.llsg.thresholding(array, threshold, mode)[source]

Array thresholding strategies.

vip_hci.psfsub.loci module

Module with a frame differencing algorithm for ADI and ADI+mSDI post-processing.

[PUE12]
Pueyo et al. 2012
Application of a Damped Locally Optimized Combination of Images Method to the Spectral Characterization of Faint Companions Using an Integral Field Spectrograph
The Astrophysical Journal Supplements, Volume 199, p. 6
class vip_hci.psfsub.loci.XLOCI_Params(cube: ndarray | None = None, angle_list: ndarray | None = None, scale_list: ndarray | None = None, fwhm: float = 4, metric: Enum = Metric.MANHATTAN, dist_threshold: int = 100, delta_rot: float | Tuple[float] = (0.1, 1), delta_sep: float | Tuple[float] = (0.1, 1), radius_int: int = 0, asize: int = 4, n_segments: int = 4, nproc: int = 1, solver: Enum = Solver.LSTSQ, tol: float = 0.01, optim_scale_fact: float = 2, adimsdi: Enum = Adimsdi.SKIPADI, imlib: Enum = Imlib.VIPFFT, interpolation: Enum = Interpolation.LANCZOS4, collapse: Enum = Collapse.MEDIAN, verbose: bool = True, full_output: bool = False)[source]

Bases: object

Set of parameters for the LOCI algorithm.

See function xloci below for the documentation.

adimsdi: Enum = 'skipadi'
angle_list: ndarray = None
asize: int = 4
collapse: Enum = 'median'
cube: ndarray = None
delta_rot: float | Tuple[float] = (0.1, 1)
delta_sep: float | Tuple[float] = (0.1, 1)
dist_threshold: int = 100
full_output: bool = False
fwhm: float = 4
imlib: Enum = 'vip-fft'
interpolation: Enum = 'lanczos4'
metric: Enum = 'manhattan'
n_segments: int = 4
nproc: int = 1
optim_scale_fact: float = 2
radius_int: int = 0
scale_list: ndarray = None
solver: Enum = 'lstsq'
tol: float = 0.01
verbose: bool = True
vip_hci.psfsub.loci.xloci(*all_args: List, **all_kwargs: dict)[source]

Locally Optimized Combination of Images (LOCI) algorithm as in [LAF07]. The PSF is modeled (for ADI and ADI+mSDI) with a least-square combination of neighbouring frames (solving the equation a x = b by computing a vector x of coefficients that minimizes the Euclidean 2-norm || b - a x ||^2).

This algorithm is also compatible with IFS data to perform LOCI-SDI, in a similar fashion as suggested in [PUE12] (albeit without dampening zones).

Parameters:
  • all_args (list, optional) – Positionnal arguments for the LOCI algorithm. Full list of parameters below.

  • all_kwargs (dictionary, optional) – Mix of keyword arguments that can initialize a LOCIParams and the optional ‘rot_options’ dictionnary, with keyword values for “border_mode”, “mask_val”, “edge_blend”, “interp_zeros”, “ker” (see documentation of vip_hci.preproc.frame_rotate). Can also contain a LOCIParams named as algo_params.

  • parameters (LOCI)

  • ----------

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

  • angle_list (numpy ndarray, 1d) – Corresponding parallactic angle for each frame.

  • scale_list (numpy ndarray, 1d, optional) – If provided, triggers mSDI reduction. These should be the scaling factors used to re-scale the spectral channels and align the speckles in case of IFS data (ADI+mSDI cube). Usually, these can be approximated by the last channel wavelength divided by the other wavelengths in the cube (more thorough approaches can be used to get the scaling factors, e.g. with vip_hci.preproc.find_scal_vector).

  • fwhm (float, optional) – Size of the FWHM in pixels. Default is 4.

  • metric (Enum, see vip_hci.config.paramenum.Metric) – Distance metric to be used (‘cityblock’, ‘cosine’, ‘euclidean’, ‘l1’, ‘l2’, ‘manhattan’, ‘correlation’, etc). It uses the scikit-learn function sklearn.metrics.pairwise.pairwise_distances (check its documentation).

  • dist_threshold (int, optional) – Indices with a distance larger than dist_threshold percentile will initially discarded. 100 by default.

  • delta_rot (float or tuple of floats, optional) – Factor for adjusting the parallactic angle threshold, expressed in FWHM. Default is 1 (excludes 1 FWHM on each side of the considered frame). If a tuple of two floats is provided, they are used as the lower and upper intervals for the threshold (grows linearly as a function of the separation).

  • 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).

  • radius_int (int, optional) – The radius of the innermost annulus. By default is 0, if >0 then the central circular region is discarded.

  • asize (int, optional) – The size of the annuli, in pixels.

  • n_segments (int or list of int or 'auto', optional) – The number of segments for each annulus. When a single integer is given it is used for all annuli. When set to ‘auto’, the number of segments is automatically determined for every annulus, based on the annulus width.

  • nproc (None or int, optional) – Number of processes for parallel computing. If None the number of processes will be set to cpu_count()/2. By default the algorithm works in single-process mode.

  • solver (Enum, see vip_hci.config.paramenum.Solver) – Choosing the solver of the least squares problem. lstsq uses the standard scipy least squares solver. nnls uses the scipy non-negative least-squares solver.

  • tol (float, optional) – Valid when solver is set to lstsq. Sets the cutoff for ‘small’ singular values; used to determine effective rank of a. Singular values smaller than tol * largest_singular_value are considered zero. Smaller values of tol lead to smaller residuals (more aggressive subtraction).

  • optim_scale_fact (float, optional) – If >1, the least-squares optimization is performed on a larger segment, similar to LOCI. The optimization segments share the same inner radius, mean angular position and angular width as their corresponding subtraction segments.

  • adimsdi (Enum, see vip_hci.config.paramenum.Adimsdi) –

    Changes the way the 4d cubes (ADI+mSDI) are processed.

    skipadi: the multi-spectral frames are rescaled wrt the largest wavelength to align the speckles and the least-squares model is subtracted on each spectral cube separately.

    double: a first subtraction is done on the rescaled spectral frames (as in the skipadi case). Then the residuals are processed again in an ADI fashion.

  • imlib (Enum, see vip_hci.config.paramenum.Imlib) – See the documentation of the vip_hci.preproc.frame_rotate function.

  • interpolation (Enum, see vip_hci.config.paramenum.Interpolation) – See the documentation of the vip_hci.preproc.frame_rotate function.

  • collapse (Enum, see vip_hci.config.paramenum.Collapse) – Sets the way of collapsing the frames for producing a final image.

  • verbose (bool, optional) – If True prints info to stdout.

  • full_output (bool, optional) – Whether to return the final median combined image only or along with 2 other residual cubes (before and after derotation).

Returns:

  • cube_res (numpy ndarray, 3d) – [full_output=True] Cube of residuals.

  • cube_der (numpy ndarray, 3d) – [full_output=True] Derotated cube of residuals.

  • frame_der_median (numpy ndarray, 2d) – Median combination of the de-rotated cube of residuals.

vip_hci.psfsub.medsub module

Implementation of a median subtraction algorithm for model PSF subtraction in high-contrast imaging sequences. Median-ADI was originally proposed in [MAR06], while median-SDI (also referred to as spectral deconvolution) was proposed in [SPA02] and further developed in [THA07].

[MAR06] (1,2,3,4)
Marois et al. 2006
Angular Differential Imaging: A Powerful High-Contrast Imaging Technique
The Astrophysical Journal, Volume 641, Issue 1, pp. 556-564
[SPA02] (1,2,3)
Sparks & Ford 2002
Imaging Spectroscopy for Extrasolar Planet Detection
The Astrophysical Journal, Volume 578, Issue 1, pp. 543-564
[THA07] (1,2,3)
Thatte et al. 2007
Very high contrast integral field spectroscopy of AB Doradus C: 9-mag contrast at 0.2arcsec without a coronagraph using spectral deconvolution
MNRAS, Volume 378, Issue 4, pp. 1229-1236
class vip_hci.psfsub.medsub.MEDIAN_SUB_Params(cube: ndarray | None = None, angle_list: ndarray | None = None, scale_list: ndarray | None = None, flux_sc_list: ndarray | None = None, fwhm: float = 4, radius_int: int = 0, asize: int = 4, delta_rot: int = 1, delta_sep: float | Tuple[float] = (0.1, 1), mode: str = 'fullfr', nframes: int = 4, sdi_only: bool = False, imlib: Enum = Imlib.VIPFFT, interpolation: Enum = Interpolation.LANCZOS4, collapse: Enum = Collapse.MEDIAN, nproc: int = 1, full_output: bool = False, verbose: bool = True)[source]

Bases: object

Set of parameters for the median subtraction module.

See function median_sub for documentation.

angle_list: ndarray = None
asize: int = 4
collapse: Enum = 'median'
cube: ndarray = None
delta_rot: int = 1
delta_sep: float | Tuple[float] = (0.1, 1)
flux_sc_list: ndarray = None
full_output: bool = False
fwhm: float = 4
imlib: Enum = 'vip-fft'
interpolation: Enum = 'lanczos4'
mode: str = 'fullfr'
nframes: int = 4
nproc: int = 1
radius_int: int = 0
scale_list: ndarray = None
sdi_only: bool = False
verbose: bool = True
vip_hci.psfsub.medsub.median_sub(*all_args: List, **all_kwargs: dict)[source]

Perform (smart) median-ADI or median-SDI.

In the case of angular differential imaging (ADI), the algorithm is based on [MAR06]. The ADI+IFS method is an extension of this basic idea to multi-spectral cubes, combining ADI with spectral deconvolution (also called spectral differential imaging or SDI).

References: [MAR06] for median-ADI; [SPA02] and [THA07] for SDI.

Parameters:
  • all_args (list, optional) – Positionnal arguments for the median_sub algorithm. Full list of parameters is provided below.

  • all_kwargs (dictionary, optional) – Mix of keyword arguments that can initialize a MEDIAN_SUB_Params and the optional rot_options dictionary (with keywords border_mode, mask_val, edge_blend, interp_zeros, ker; see docstrings of vip_hci.preproc.frame_rotate). Can also contain a MEDIAN_SUB_Params object/dictionary named algo_params.

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

  • angle_list (numpy ndarray, 1d) – Corresponding parallactic angle for each frame.

  • scale_list (numpy ndarray, 1d, optional) – If provided, triggers mSDI reduction. These should be the scaling factors used to re-scale the spectral channels and align the speckles in case of IFS data (ADI+mSDI cube). Usually, these can be approximated by the last channel wavelength divided by the other wavelengths in the cube (more thorough approaches can be used to get the scaling factors, e.g. with vip_hci.preproc.find_scal_vector).

  • flux_sc_list (numpy ndarray, 1d) – In the case of IFS data (ADI+SDI), this is the list of flux scaling factors applied to each spectral frame after geometrical rescaling. These should be set to either the ratio of stellar fluxes between the last spectral channel and the other channels, or to the second output of preproc.find_scal_vector (when using 2 free parameters). If not provided, the algorithm will still work, but with a lower efficiency at subtracting the stellar halo.

  • fwhm (float or 1d numpy array) – Known size of the FWHM in pixels to be used. Default is 4.

  • 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 pixels.

  • delta_rot (int, optional) – Factor for increasing the parallactic angle threshold, expressed in FWHM. Default is 1 (excludes 1 FWHM 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).

  • mode ({'fullfr', 'annular'}, str optional) – In fullfr mode only the median frame is subtracted, in annular mode also the 4 closest frames given a PA threshold (annulus-wise) are subtracted.

  • nframes (int or None, optional) – Number of frames (even value) to be used for building the optimized reference PSF when working in annular mode. None by default, which means that all frames, excluding the thresholded ones, are used.

  • sdi_only (bool, optional) – In the case of IFS data (ADI+SDI), whether to perform median-SDI, or median-ASDI (default).

  • imlib (Enum, see vip_hci.config.paramenum.Imlib) – See the documentation of vip_hci.preproc.frame_rotate.

  • interpolation (Enum, see vip_hci.config.paramenum.Interpolation) – See the documentation of the vip_hci.preproc.frame_rotate function.

  • collapse (Enum, see vip_hci.config.paramenum.Collapse) – Sets how temporal residual frames should be combined to produce an ADI image.

  • nproc (None or int, optional) – Number of processes for parallel computing. If None the number of processes will be set to cpu_count()/2. By default the algorithm works in single-process mode.

  • full_output (bool, optional) – Whether to return the final median combined image only or with other intermediate arrays.

  • verbose (bool, optional) – If True prints to stdout intermediate info.

Returns:

  • cube_out (numpy ndarray, 3d) – [full_output=True] The cube of residuals.

  • cube_der (numpy ndarray, 3d) – [full_output=True] The derotated cube of residuals.

  • frame (numpy ndarray, 2d) – Median combination of the de-rotated cube.

vip_hci.psfsub.nmf_fullfr module

Module with PSF reference approximation using Non-negative matrix factorization for ADI and RDI data, in full frames.

[LEE99] (1,2)
Lee & Seung 1999
Learning the parts of objects by non-negative matrix factorization
Nature, Volume 401, Issue 6755, pp. 788-791
class vip_hci.psfsub.nmf_fullfr.NMF_Params(cube: ~numpy.ndarray | None = None, angle_list: ~numpy.ndarray | None = None, cube_ref: ~numpy.ndarray | None = None, ncomp: int = 1, scaling: ~enum.Enum | None = None, max_iter: int = 10000, random_state: int | None = None, mask_center_px: int | None = None, source_xy: ~typing.Tuple[int] | None = None, delta_rot: float = 1, fwhm: float = 4, init_svd: ~enum.Enum = Initsvd.NNDSVD, collapse: ~enum.Enum = Collapse.MEDIAN, full_output: bool = False, verbose: bool = True, cube_sig: ~numpy.ndarray | None = None, handle_neg: ~enum.Enum = HandleNeg.MASK, nmf_args: dict = <factory>)[source]

Bases: object

Set of parameters for the NMF full-frame algorithm.

See function nmf below for the documentation.

angle_list: ndarray = None
collapse: Enum = 'median'
cube: ndarray = None
cube_ref: ndarray = None
cube_sig: ndarray = None
delta_rot: float = 1
full_output: bool = False
fwhm: float = 4
handle_neg: Enum = 'mask'
init_svd: Enum = 'nndsvd'
mask_center_px: int = None
max_iter: int = 10000
ncomp: int = 1
nmf_args: dict
random_state: int = None
scaling: Enum = None
source_xy: Tuple[int] = None
verbose: bool = True
vip_hci.psfsub.nmf_fullfr.nmf(*all_args: List, **all_kwargs: dict)[source]

Non Negative Matrix Factorization [LEE99] for ADI sequences [GOM17]. Alternative to the full-frame ADI-PCA processing that does not rely on SVD or ED for obtaining a low-rank approximation of the datacube. This function embeds the scikit-learn NMF algorithm solved through either the coordinate descent or the multiplicative update method.

Parameters:
  • all_args (list, optional) – Positionnal arguments for the NMF algorithm. Full list of parameters below.

  • all_kwargs (dictionary, optional) – Mix of keyword arguments that can initialize a NMFParams and the optional ‘rot_options’ dictionnary, with keyword values for “border_mode”, “mask_val”, “edge_blend”, “interp_zeros”, “ker” (see documentation of vip_hci.preproc.frame_rotate). Can also contain a NMFParams named as algo_params.

  • parameters (NMF)

  • ----------

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

  • angle_list (numpy ndarray, 1d) – Corresponding parallactic angle for each frame.

  • cube_ref (numpy ndarray, 3d, optional) – Reference library cube. For Reference Star Differential Imaging.

  • ncomp (int, optional) – How many components are used as for low-rank approximation of the datacube.

  • scaling (Enum, see vip_hci.config.paramenum.Scaling) – 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.

  • max_iter (int optional) – The number of iterations for the coordinate descent solver.

  • random_state (int or None, optional) – Controls the seed for the Pseudo Random Number generator.

  • 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, with source_xy as the coordinates X,Y of the center of the annulus where the PA criterion is estimated. When ncomp is a tuple, a PCA grid is computed and the S/Ns (mean value in a 1xFWHM circular aperture) of the given (X,Y) coordinates are computed.

  • delta_rot (float, optional) – Factor for tunning the parallactic angle threshold, expressed in FWHM. Default is 1 (excludes 1xFWHM on each side of the considered frame).

  • fwhm (float, optional) – Known size of the FWHM in pixels to be used. Default value is 4.

  • init_svd (Enum, see vip_hci.config.paramenum.Initsvd) – Method used to initialize the iterative procedure to find H and W. ‘nndsvd’: non-negative double SVD recommended for sparseness ‘nndsvda’: NNDSVD where zeros are filled with the average of cube; recommended when sparsity is not desired ‘random’: random initial non-negative matrix

  • collapse (Enum, see vip_hci.config.paramenum.Collapse) – Sets the way of collapsing the frames for producing a final image.

  • full_output (boolean, optional) – Whether to return the final median combined image only or with other intermediate arrays.

  • verbose ({True, False}, bool optional) – If True prints intermediate info and timing.

  • handle_neg (Enum, see vip_hci.config.paramenum.HandleNeg) – Determines how to handle negative values: mask them, set them to zero, or subtract the minimum value in the arrays. Note: ‘mask’ or ‘null’ may leave significant artefacts after derotation of residual cube => those options should be used carefully (e.g. with proper treatment of masked values in non-derotated cube of residuals).

  • nmf_args (dictionary, optional) – Additional arguments for scikit-learn NMF algorithm. See: https://scikit-learn.org/stable/modules/generated/sklearn.decomposition.NMF.html

Returns:

  • If full_output is False the final frame is returned. If True the algorithm

  • returns the reshaped NMF components, the reconstructed cube, the residuals,

  • the residuals derotated and the final frame.

vip_hci.psfsub.nmf_local module

Module with NMF algorithm in concentric annuli for ADI/RDI.

class vip_hci.psfsub.nmf_local.NMF_ANNULAR_Params(cube: ~numpy.ndarray | None = None, angle_list: ~numpy.ndarray | None = None, cube_ref: ~numpy.ndarray | None = None, radius_int: int = 0, fwhm: float = 4, asize: int = 4, n_segments: int = 1, delta_rot: float | ~typing.Tuple[float] = (0.1, 1), ncomp: int = 1, init_svd: ~enum.Enum = Initsvd.NNDSVD, nproc: int = 1, min_frames_lib: int = 2, max_frames_lib: int = 200, scaling: ~enum.Enum | None = None, imlib: ~enum.Enum = Imlib.VIPFFT, interpolation: ~enum.Enum = Interpolation.LANCZOS4, collapse: ~enum.Enum = Collapse.MEDIAN, full_output: bool = False, verbose: bool = True, theta_init: float = 0, weights: ~typing.List | None = None, cube_sig: ~numpy.ndarray | None = None, handle_neg: ~enum.Enum = HandleNeg.MASK, max_iter: int = 1000, random_state: int | None = None, nmf_args: dict = <factory>)[source]

Bases: object

Set of parameters for the NMF annular algorithm.

See function nmf_annular below for the documentation.

angle_list: ndarray = None
asize: int = 4
collapse: Enum = 'median'
cube: ndarray = None
cube_ref: ndarray = None
cube_sig: ndarray = None
delta_rot: float | Tuple[float] = (0.1, 1)
full_output: bool = False
fwhm: float = 4
handle_neg: Enum = 'mask'
imlib: Enum = 'vip-fft'
init_svd: Enum = 'nndsvd'
interpolation: Enum = 'lanczos4'
max_frames_lib: int = 200
max_iter: int = 1000
min_frames_lib: int = 2
n_segments: int = 1
ncomp: int = 1
nmf_args: dict
nproc: int = 1
radius_int: int = 0
random_state: int = None
scaling: Enum = None
theta_init: float = 0
verbose: bool = True
weights: List = None
vip_hci.psfsub.nmf_local.nmf_annular(*all_args: List, **all_kwargs: dict)[source]

Non Negative Matrix Factorization in concentric annuli, for ADI/RDI sequences. Alternative to the annular ADI-PCA processing that does not rely on SVD or ED for obtaining a low-rank approximation of the datacube. This function embeds the scikit-learn NMF algorithm solved through either the coordinate descent or the multiplicative update method.

Parameters:
  • all_args (list, optional) – Positionnal arguments for the NMF annular algorithm. Full list of parameters below.

  • all_kwargs (dictionary, optional) – Mix of keyword arguments that can initialize a NMFAnnParams and the optional ‘rot_options’ dictionnary, with keyword values for “border_mode”, “mask_val”, “edge_blend”, “interp_zeros”, “ker” (see documentation of vip_hci.preproc.frame_rotate). Can also contain a NMFAnnParams named as algo_params..

  • parameters (NMF annular)

  • ----------

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

  • angle_list (numpy ndarray, 1d) – Corresponding parallactic angle for each frame.

  • cube_ref (numpy ndarray, 3d, optional) – Reference library cube. For Reference Star Differential Imaging.

  • radius_int (int, optional) – The radius of the innermost annulus. By default is 0, if >0 then the central circular region is discarded.

  • fwhm (float, optional) – Size of the FWHM in pixels. Default is 4.

  • asize (float, optional) – The size of the annuli, in pixels.

  • n_segments (int or list of ints or 'auto', optional) – The number of segments for each annulus. When a single integer is given it is used for all annuli. When set to ‘auto’, the number of segments is automatically determined for every annulus, based on the annulus width.

  • delta_rot (int, optional) – Factor for adjusting the parallactic angle threshold, expressed in FWHM. Default is 1 (excludes 1 FWHM on each side of the considered frame). If a tuple of two floats is provided, they are used as the lower and upper intervals for the threshold (grows linearly as a function of the separation). !!! Important: this is used even if a reference cube is provided for RDI. This is to allow ARDI (PCA library built from both science and reference cubes). If you want to do pure RDI, set delta_rot to an arbitrarily high value such that the condition is never fulfilled for science frames to make it in the PCA library.

  • ncomp (int, optional) – How many components are used as for low-rank approximation of the datacube.

  • scaling ({None, "temp-mean", spat-mean", "temp-standard",) –

    “spat-standard”}, None or str optional Pixel-wise scaling mode using sklearn.preprocessing.scale function. If set to None, the input matrix is left untouched. Otherwise:

    • temp-mean: temporal px-wise mean is subtracted.

    • spat-mean: spatial mean is subtracted.

    • temp-standard: temporal mean centering plus scaling pixel values to unit variance (temporally).

    • spat-standard: spatial mean centering plus scaling pixel values to unit variance (spatially).

    DISCLAIMER: Using temp-mean or temp-standard scaling can improve the speckle subtraction for ASDI or (A)RDI reductions. Nonetheless, this involves a sort of c-ADI preprocessing, which (i) can be dangerous for datasets with low amount of rotation (strong self-subtraction), and (ii) should probably be referred to as ARDI (i.e. not RDI stricto sensu).

  • max_iter (int optional) – The number of iterations for the coordinate descent solver.

  • random_state (int or None, optional) – Controls the seed for the Pseudo Random Number generator.

  • 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, with source_xy as the coordinates X,Y of the center of the annulus where the PA criterion is estimated. When ncomp is a tuple, a PCA grid is computed and the S/Ns (mean value in a 1xFWHM circular aperture) of the given (X,Y) coordinates are computed.

  • delta_rot – Factor for tunning the parallactic angle threshold, expressed in FWHM. Default is 1 (excludes 1xFWHM on each side of the considered frame).

  • init_svd (str, optional {'nnsvd','nnsvda','random'}) – Method used to initialize the iterative procedure to find H and W. ‘nndsvd’: non-negative double SVD recommended for sparseness ‘nndsvda’: NNDSVD where zeros are filled with the average of cube; recommended when sparsity is not desired ‘random’: random initial non-negative matrix

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

  • full_output (boolean, optional) – Whether to return the final median combined image only or with other intermediate arrays.

  • verbose ({True, False}, bool optional) – If True prints intermediate info and timing.

  • nmf_args (dictionary, optional) – Additional arguments for scikit-learn NMF algorithm. See: https://scikit-learn.org/stable/modules/generated/sklearn.decomposition.NMF.html

Returns:

  • If full_output is False the final frame is returned. If True the algorithm

  • returns the reshaped NMF components, the reconstructed cube, the residuals,

  • the residuals derotated and the final frame.

vip_hci.psfsub.pca_fullfr module

Full-frame PCA algorithm for ADI, (ADI+)RDI and (ADI+)mSDI (IFS data) cubes.

Options :

  • Full-frame PCA, using the whole cube as the PCA reference library in the case of ADI or ADI+mSDI (IFS cube), or a sequence of reference frames (reference star) in the case of RDI. For ADI a big data matrix NxP, where N is the number of frames and P the number of pixels in a frame is created. Then PCA is done through eigen-decomposition of the covariance matrix (~$DD^T$) or the SVD of the data matrix. SVD can be calculated using different libraries including the fast randomized SVD.

  • Full-frame incremental PCA for big (larger than available memory) cubes.

[AMA12] (1,2,3)
Amara & Quanz 2012
PYNPOINT: an image processing package for finding exoplanets
MNRAS, Volume 427, Issue 1, pp. 948-955
[CHR19] (1,2)
Christiaens et al. 2019
Separating extended disc features from the protoplanet in PDS 70 using VLT/SINFONI
MNRAS, Volume 486, Issue 4, pp. 5819-5837
[HAL09] (1,2)
Halko et al. 2009
Finding structure with randomness: Probabilistic algorithms for constructing approximate matrix decompositions
arXiv e-prints
[REN23]
Ren 2023
Karhunen-Loève data imputation in high-contrast imaging
Astronomy & Astrophysics, Volume 679, p. 8
class vip_hci.psfsub.pca_fullfr.PCA_Params(cube: ndarray | None = None, angle_list: ndarray | None = None, cube_ref: ndarray | None = None, scale_list: ndarray | None = None, ncomp: Tuple | List | float | int = 1, svd_mode: Enum = SvdMode.LAPACK, scaling: Enum | None = None, mask_center_px: int | None = None, source_xy: Tuple[int] | None = None, delta_rot: int | None = None, fwhm: float = 4, adimsdi: Enum = Adimsdi.SINGLE, crop_ifs: bool = True, imlib: Enum = Imlib.VIPFFT, imlib2: Enum = Imlib.VIPFFT, interpolation: Enum = Interpolation.LANCZOS4, collapse: Enum = Collapse.MEDIAN, collapse_ifs: Enum = Collapse.MEAN, ifs_collapse_range: str | Tuple[int] = 'all', mask_rdi: ndarray | None = None, check_memory: bool = True, batch: int | float | None = None, nproc: int = 1, full_output: bool = False, verbose: bool = True, weights: ndarray | None = None, left_eigv: bool = False, min_frames_pca: int = 10, cube_sig: ndarray | None = None)[source]

Bases: object

Set of parameters for the PCA module.

See function pca below for the documentation.

adimsdi: Enum = 'single'
angle_list: ndarray = None
batch: int | float = None
check_memory: bool = True
collapse: Enum = 'median'
collapse_ifs: Enum = 'mean'
crop_ifs: bool = True
cube: ndarray = None
cube_ref: ndarray = None
cube_sig: ndarray = None
delta_rot: int = None
full_output: bool = False
fwhm: float = 4
ifs_collapse_range: str | Tuple[int] = 'all'
imlib: Enum = 'vip-fft'
imlib2: Enum = 'vip-fft'
interpolation: Enum = 'lanczos4'
left_eigv: bool = False
mask_center_px: int = None
mask_rdi: ndarray = None
min_frames_pca: int = 10
ncomp: Tuple | List | float | int = 1
nproc: int = 1
scale_list: ndarray = None
scaling: Enum = None
source_xy: Tuple[int] = None
svd_mode: Enum = 'lapack'
verbose: bool = True
weights: ndarray = None
vip_hci.psfsub.pca_fullfr.pca(*all_args: List, **all_kwargs: dict)[source]

Full-frame PCA algorithm applied to PSF substraction.

The reference PSF and the quasi-static speckle pattern are modeled using Principal Component Analysis. Depending on the input parameters this PCA function can work in ADI, RDI or mSDI (IFS data) mode.

ADI: the target cube itself is used to learn the PCs and to obtain a low-rank approximation model PSF (star + speckles). Both cube_ref` and scale_list must be None. The full-frame ADI-PCA implementation is based on [AMA12] and [SOU12]. If batch is provided then the cube is processed with incremental PCA as described in [GOM17].

(ADI+)RDI: if a reference cube is provided (cube_ref), its PCs are used to reconstruct the target frames to obtain the model PSF (star + speckles).

(ADI+)mSDI (IFS data): if a scaling vector is provided (scale_list) and the cube is a 4d array [# channels, # adi-frames, Y, X], it’s assumed it contains several multi-spectral frames acquired in pupil-stabilized mode. A single or two stages PCA can be performed, depending on adimsdi, as explained in [CHR19].

Parameters:
  • all_args (list, optional) – Positionnal arguments for the PCA algorithm. Full list of parameters below.

  • all_kwargs (dictionary, optional) – Mix of keyword arguments that can initialize a PCA_Params and the optional ‘rot_options’ dictionnary (with keyword values border_mode, mask_val, edge_blend, interp_zeros, ker; see docstring of vip_hci.preproc.frame_rotate). Can also contain a PCA_Params dictionary named algo_params.

  • parameters (PCA)

  • ----------

  • cube (str or numpy ndarray, 3d or 4d) – Input cube (ADI or ADI+mSDI). If 4D, the first dimension should be spectral. If a string is given, it must correspond to the path to the fits file to be opened in memmap mode (incremental PCA-ADI of 3D cubes only).

  • angle_list (numpy ndarray, 1d) – Corresponding parallactic angle for each frame.

  • cube_ref (3d or 4d numpy ndarray, or list of 3D numpy ndarray, optional) – Reference library cube for Reference Star Differential Imaging. Should be 3D, except if input cube is 4D and no scale_list is provided, reference cube can then either be 4D or a list of 3D cubes (i.e. providing the reference cube for each individual spectral cube).

  • scale_list (numpy ndarray, 1d, optional) – If provided, triggers mSDI reduction. These should be the scaling factors used to re-scale the spectral channels and align the speckles in case of IFS data (ADI+mSDI cube). Usually, the scaling factors are the last channel wavelength divided by the other wavelengths in the cube (more thorough approaches can be used to get the scaling factors, e.g. with vip_hci.preproc.find_scal_vector).

  • ncomp (int, float, tuple of int/None, or list, optional) –

    How many PCs are used as a lower-dimensional subspace to project the target frames.

    • ADI (cube is a 3d array): if an int is provided, ncomp is the

    number of PCs extracted from cube itself. If ncomp is a float in the interval [0, 1] then it corresponds to the desired cumulative explained variance ratio (the corresponding number of components is estimated). If ncomp is a tuple of two integers, then it corresponds to an interval of PCs in which final residual frames are computed (optionally, if a tuple of 3 integers is passed, the third value is the step). If ncomp is a list of int, these will be used to calculate residual frames. When ncomp is a tuple or list, and source_xy is not None, then the S/Ns (mean value in a 1xFWHM circular aperture) of the given (X,Y) coordinates are computed.

    • ADI+RDI (cube and cube_ref are 3d arrays): ncomp is the

    number of PCs obtained from cube_ref. If ncomp is a tuple, then it corresponds to an interval of PCs (obtained from cube_ref) in which final residual frames are computed. If ncomp is a list of int, these will be used to calculate residual frames. When ncomp is a tuple or list, and source_xy is not None, then the S/Ns (mean value in a 1xFWHM circular aperture) of the given (X,Y) coordinates are computed.

    • ADI or ADI+RDI (cube is a 4d array): same input format allowed as

    above. If ncomp is a list with the same length as the number of channels, each element of the list will be used as ncomp value (be it int, float or tuple) for each spectral channel. If not a list or a list with a different length as the number of spectral channels, these will be tested for all spectral channels respectively.

    • ADI+mSDI (cube is a 4d array and adimsdi="single"): ncomp

    is the number of PCs obtained from the whole set of frames (n_channels * n_adiframes). If ncomp is a float in the interval (0, 1] then it corresponds to the desired CEVR, and the corresponding number of components will be estimated. If ncomp is a tuple, then it corresponds to an interval of PCs in which final residual frames are computed. If ncomp is a list of int, these will be used to calculate residual frames. When ncomp is a tuple or list, and source_xy is not None, then the S/Ns (mean value in a 1xFWHM circular aperture) of the given (X,Y) coordinates are computed.

    • ADI+mSDI (cube is a 4d array and adimsdi="double"): ncomp

    must be a tuple, where the first value is the number of PCs obtained from each multi-spectral frame (if None then this stage will be skipped and the spectral channels will be combined without subtraction); the second value sets the number of PCs used in the second PCA stage, ADI-like 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 (Enum, see vip_hci.config.paramenum.SvdMode) – Switch for the SVD method/library to be used.

  • scaling (Enum, see vip_hci.config.paramenum.Scaling) – Pixel-wise scaling mode using sklearn.preprocessing.scale function. If set to None, the input matrix is left untouched.

  • 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, with source_xy as the coordinates X,Y of the center of the annulus where the PA criterion is estimated. When ncomp is a tuple, a PCA grid is computed and the S/Ns (mean value in a 1xFWHM circular aperture) of the given (X,Y) coordinates are computed.

  • delta_rot (int, optional) – Factor for tuning the parallactic angle threshold, expressed in FWHM. Default is 1 (excludes 1xFWHM on each side of the considered frame).

  • fwhm (float, list or 1d numpy array, optional) – Known size of the FWHM in pixels to be used. Default value is 4. Can be a list or 1d numpy array for a 4d input cube with no scale_list.

  • adimsdi (Enum, see vip_hci.config.paramenum.Adimsdi) – Changes the way the 4d cubes (ADI+mSDI) are processed. Basically it determines whether a single or double pass PCA is going to be computed.

  • crop_ifs (bool, optional) – [adimsdi=’single’] If True cube is cropped at the moment of frame rescaling in wavelength. This is recommended for large FOVs such as the one of SPHERE, but can remove significant amount of information close to the edge of small FOVs (e.g. SINFONI).

  • imlib (Enum, see vip_hci.config.paramenum.Imlib) – See the documentation of vip_hci.preproc.frame_rotate.

  • imlib2 (Enum, see vip_hci.config.paramenum.Imlib) – See the documentation of vip_hci.preproc.cube_rescaling_wavelengths.

  • interpolation (Enum, see vip_hci.config.paramenum.Interpolation) – See the documentation of the vip_hci.preproc.frame_rotate function.

  • collapse (Enum, see vip_hci.config.paramenum.Collapse) – Sets how temporal residual frames should be combined to produce an ADI image.

  • collapse_ifs (Enum, see vip_hci.config.paramenum.Collapse) – Sets how spectral residual frames should be combined to produce an mSDI image.

  • ifs_collapse_range (str 'all' or tuple of 2 int) – If a tuple, it should contain the first and last channels where the mSDI residual channels will be collapsed (by default collapses all channels).

  • mask_rdi (tuple of two numpy array or one signle 2d numpy array, opt) – If provided, binary mask(s) will be used either in RDI mode or in ADI+mSDI (2 steps) mode. They will be used as boat and anchor masks following the procedure described in [REN23], which is useful to avoid self-subtraction in the presence of a bright disc signal. If only one mask is provided, the boat images will be unmasked (i.e., full frames will be used).

  • check_memory (bool, optional) – If True, it checks that the input cube is smaller than the available system memory.

  • batch (None, int or float, optional) – When it is not None, it triggers the incremental PCA (for ADI and ADI+mSDI cubes). If an int is given, it corresponds to the number of frames in each sequential mini-batch. If a float (0, 1] is given, it corresponds to the size of the batch is computed wrt the available memory in the system.

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

  • full_output (bool, optional) – Whether to return the final median combined image only or with other intermediate arrays.

  • verbose (bool, optional) – If True prints intermediate info and timing.

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

  • left_eigv (bool, optional) – Whether to use rather left or right singularvectors This mode is not compatible with ‘mask_rdi’ and ‘batch’

  • min_frames_pca (int, optional) – Minimum number of frames required in the PCA library. An error is raised if less than such number of frames can be found.

  • cube_sig (numpy ndarray, opt) – Cube with estimate of significant authentic signals. If provided, this will be subtracted before projecting considering the science cube as reference cube.

Returns:

  • final_residuals_cube (List of numpy ndarray) – [(ncomp is tuple or list) & (source_xy=None or full_output=True)] List of residual final PCA frames obtained for a grid of PC values.

  • frame (numpy ndarray) – [ncomp is scalar or source_xy!=None] 2D array, median combination of the de-rotated/re-scaled residuals cube.

  • pcs (numpy ndarray) – [full_output=True, source_xy=None] Principal components. Valid for ADI cubes 3D or 4D (i.e. scale_list=None). This is also returned when batch is not None (incremental PCA).

  • recon_cube, recon (numpy ndarray) – [full_output=True] Reconstructed cube. Valid for ADI cubes 3D or 4D (i.e. scale_list=None)

  • residuals_cube (numpy ndarray) – [full_output=True] Residuals cube. Valid for ADI cubes 3D or 4D (i.e. scale_list=None)

  • residuals_cube_ (numpy ndarray) – [full_output=True] Derotated residuals cube. Valid for ADI cubes 3D or 4D (i.e. scale_list=None)

  • residuals_cube_channels (numpy ndarray) – [full_output=True, adimsdi=’double’] Residuals for each multi-spectral cube. Valid for ADI+mSDI (4D) cubes (when scale_list is provided)

  • residuals_cube_channels_ (numpy ndarray) – [full_output=True, adimsdi=’double’] Derotated final residuals. Valid for ADI+mSDI (4D) cubes (when scale_list is provided)

  • cube_allfr_residuals (numpy ndarray) – [full_output=True, adimsdi=’single’] Residuals cube (of the big cube with channels and time processed together). Valid for ADI+mSDI (4D) cubes (when scale_list is provided)

  • cube_adi_residuals (numpy ndarray) – [full_output=True, adimsdi=’single’] Residuals cube (of the big cube with channels and time processed together) after de-scaling the wls. Valid for ADI+mSDI (4D) (when scale_list is provided).

  • medians (numpy ndarray) – [full_output=True, source_xy=None] This is also returned when batch is not None (incremental PCA).

  • ifs_adi_frames (numpy ndarray) – [full_output=True, 4D input cube, scale_list=None] This is the cube of individual ADI reductions for each channel of the IFS cube.

vip_hci.psfsub.pca_local module

Module with local/smart PCA (annulus or patch-wise in a multi-processing fashion) model PSF subtraction for ADI, ADI+SDI (IFS) and ADI+RDI datasets.

[ABS13]
Absil et al. 2013
Searching for companions down to 2 AU from beta Pictoris using the L’-band AGPM coronagraph on VLT/NACO
Astronomy & Astrophysics, Volume 559, Issue 1, p. 12
class vip_hci.psfsub.pca_local.PCA_ANNULAR_Params(cube: ndarray = None, angle_list: ndarray = None, cube_ref: ndarray = None, scale_list: ndarray = None, radius_int: int = 0, fwhm: float = 4, asize: float = 4, n_segments: int | List[int] | auto = 1, delta_rot: float | Tuple[float] = (0.1, 1), delta_sep: float | Tuple[float] = (0.1, 1), ncomp: int | Tuple | ndarray | auto = 1, svd_mode: Enum = SvdMode.LAPACK, nproc: int = 1, min_frames_lib: int = 2, max_frames_lib: int = 200, tol: float = 0.1, scaling: Enum = None, imlib: Enum = Imlib.VIPFFT, interpolation: Enum = Interpolation.LANCZOS4, collapse: Enum = Collapse.MEDIAN, collapse_ifs: Enum = Collapse.MEAN, ifs_collapse_range: all | Tuple[int] = 'all', theta_init: int = 0, weights: ndarray = None, cube_sig: ndarray = None, full_output: bool = False, verbose: bool = True, left_eigv: bool = False)[source]

Bases: object

Set of parameters for the annular PCA module.

angle_list: ndarray = None
asize: float = 4
collapse: Enum = 'median'
collapse_ifs: Enum = 'mean'
cube: ndarray = None
cube_ref: ndarray = None
cube_sig: ndarray = None
delta_rot: float | Tuple[float] = (0.1, 1)
delta_sep: float | Tuple[float] = (0.1, 1)
full_output: bool = False
fwhm: float = 4
ifs_collapse_range: all | Tuple[int] = 'all'
imlib: Enum = 'vip-fft'
interpolation: Enum = 'lanczos4'
left_eigv: bool = False
max_frames_lib: int = 200
min_frames_lib: int = 2
n_segments: int | List[int] | auto = 1
ncomp: int | Tuple | ndarray | auto = 1
nproc: int = 1
radius_int: int = 0
scale_list: ndarray = None
scaling: Enum = None
svd_mode: Enum = 'lapack'
theta_init: int = 0
tol: float = 0.1
verbose: bool = True
weights: ndarray = None
vip_hci.psfsub.pca_local.pca_annular(*all_args: List, **all_kwargs: dict)[source]

PCA model PSF subtraction for ADI, ADI+RDI or ADI+mSDI (IFS) data.

The PCA model is computed locally in each annulus (or annular sectors according to n_segments). For each sector we discard reference frames taking into account a parallactic angle threshold (delta_rot) and optionally a radial movement threshold (delta_sep) for 4d cubes.

For ADI+RDI data, it computes the principal components from the reference library/cube, forcing pixel-wise temporal standardization. The number of principal components can be automatically adjusted by the algorithm by minimizing the residuals inside each patch/region.

References: [AMA12] for PCA-ADI; [ABS13] for PCA-ADI in concentric annuli considering a parallactic angle threshold; [CHR19] for PCA-ASDI and PCA-SADI in one or two steps.

Parameters:
  • all_args (list, optional) – Positionnal arguments for the PCA annular algorithm. Full list of parameters below.

  • all_kwargs (dictionary, optional) – Mix of keyword arguments that can initialize a PCA_ANNULAR_Params and the optional ‘rot_options’ dictionnary, with keyword values for “border_mode”, “mask_val”, “edge_blend”, “interp_zeros”, “ker” (see documentation of vip_hci.preproc.frame_rotate). Can also contain a PCA_ANNULAR_Params named as algo_params.

  • parameters (PCA annular)

  • ----------

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

  • angle_list (numpy ndarray, 1d) – Corresponding parallactic angle for each frame.

  • cube_ref (numpy ndarray, 3d, optional) – Reference library cube. For Reference Star Differential Imaging.

  • scale_list (numpy ndarray, 1d, optional) – If provided, triggers mSDI reduction. These should be the scaling factors used to re-scale the spectral channels and align the speckles in case of IFS data (ADI+mSDI cube). Usually, these can be approximated by the last channel wavelength divided by the other wavelengths in the cube (more thorough approaches can be used to get the scaling factors, e.g. with vip_hci.preproc.find_scal_vector).

  • radius_int (int, optional) – The radius of the innermost annulus. By default is 0, if >0 then the central circular region is discarded.

  • fwhm (float, optional) – Size of the FWHM in pixels. Default is 4.

  • asize (float, optional) – The size of the annuli, in pixels.

  • n_segments (int or list of ints or 'auto', optional) – The number of segments for each annulus. When a single integer is given it is used for all annuli. When set to ‘auto’, the number of segments is automatically determined for every annulus, based on the annulus width.

  • delta_rot (float or tuple of floats, optional) – Factor for adjusting the parallactic angle threshold, expressed in FWHM. Default is 1 (excludes 1 FWHM on each side of the considered frame). If a tuple of two floats is provided, they are used as the lower and upper intervals for the threshold (grows linearly as a function of the separation).

  • 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).

  • ncomp ('auto', int, tuple/1d numpy array of int, list, tuple of lists, opt) –

    How many PCs are used as a lower-dimensional subspace to project the target (sectors of) frames. Depends on the dimensionality of cube.

    • ADI and ADI+RDI (cube is a 3d array): if a single integer is

    provided, then the same number of PCs will be subtracted at each separation (annulus). If a tuple is provided, then a different number of PCs will be used for each annulus (starting with the innermost one). If ncomp is set to auto then the number of PCs are calculated for each region/patch automatically. If a list of int is provided, several npc will be tried at once, but the same value of npc will be used for all annuli. If a tuple of lists of int is provided, the length of tuple should match the number of annuli and different sets of npc will be calculated simultaneously for each annulus, with the exact values of npc provided in the respective lists.

    • ADI or ADI+RDI (cube is a 4d array): same input format allowed as

    above, but with a slightly different behaviour if ncomp is a list: if it has the same length as the number of channels, each element of the list will be used as ncomp value (whether int, float or tuple) for each spectral channel. Otherwise the same behaviour as above is assumed.

    • ADI+mSDI case: ncomp must be a tuple of two integers or a list of

    tuples of two integers, with the number of PCs obtained from each multi-spectral frame (for each sector) and 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 (Enum, see vip_hci.config.paramenum.SvdMode) – Switch for the SVD method/library to be used.

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

  • min_frames_lib (int, optional) – Minimum number of frames in the PCA reference library.

  • max_frames_lib (int, optional) – Maximum number of frames in the PCA reference library. The more distant/decorrelated frames are removed from the library.

  • tol (float, optional) – Stopping criterion for choosing the number of PCs when ncomp is None. Lower values will lead to smaller residuals and more PCs.

  • scaling (Enum, see vip_hci.config.paramenum.Scaling) – Pixel-wise scaling mode using sklearn.preprocessing.scale function. If set to None, the input matrix is left untouched.

  • imlib (Enum, see vip_hci.config.paramenum.Imlib) – See the documentation of the vip_hci.preproc.frame_rotate function.

  • interpolation (Enum, see vip_hci.config.paramenum.Interpolation) – See the documentation of the vip_hci.preproc.frame_rotate function.

  • collapse (Enum, see vip_hci.config.paramenum.Collapse) – Sets the way of collapsing the frames for producing a final image.

  • collapse_ifs (Enum, see vip_hci.config.paramenum.Collapse) – Sets how spectral residual frames should be combined to produce an mSDI image.

  • ifs_collapse_range (str 'all' or tuple of 2 int) – If a tuple, it should contain the first and last channels where the mSDI residual channels will be collapsed (by default collapses all channels).

  • full_output (boolean, optional) – Whether to return the final median combined image only or with other intermediate arrays.

  • verbose (bool, optional) – If True prints to stdout intermediate info.

  • theta_init (int) – Initial azimuth [degrees] of the first segment, counting from the positive x-axis counterclockwise (irrelevant if n_segments=1).

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

  • cube_sig (numpy ndarray, opt) – Cube with estimate of significant authentic signals. If provided, this will be subtracted before projecting cube onto reference cube.

Returns:

  • frame (numpy ndarray, 2d) – [full_output=False] Median combination of the de-rotated cube.

  • array_out (numpy ndarray, 3d or 4d) – [full_output=True] Cube of residuals.

  • array_der (numpy ndarray, 3d or 4d) – [full_output=True] Cube residuals after de-rotation.

  • frame (numpy ndarray, 2d) – [full_output=True] Median combination of the de-rotated cube.

vip_hci.psfsub.rollsub module

Implementation of a roll subtraction algorithm for PSF subtraction in imaging sequences obtained with space-based instruments (e.g. JWST or HST) with different roll angles. The concept was proposed in [SCH98] for application to HST/NICMOS observations.

[SCH98] (1,2,3)
Schneider et al. 1998
**Exploration of the environments of nearby stars with the NICMOS

coronagraph: instrumental performance considerations** | Proc. SPIE Vol. 3356, pp. 222-233

class vip_hci.psfsub.rollsub.ROLL_SUB_Params(cube: ndarray | None = None, angle_list: ndarray | None = None, mode: str = 'mean', imlib: Enum = Imlib.VIPFFT, interpolation: Enum = Interpolation.LANCZOS4, collapse: Enum = Collapse.MEAN, fwhm_lp_bef: float = 0.0, fwhm_lp_aft: float = 0.0, mask_rad: float = 0.0, cube_sig: ndarray | None = None, nproc: int = 1, full_output: bool = False, verbose: bool = True)[source]

Bases: object

Set of parameters for the roll subtraction module.

See function roll_sub for documentation.

angle_list: ndarray = None
collapse: Enum = 'mean'
cube: ndarray = None
cube_sig: ndarray = None
full_output: bool = False
fwhm_lp_aft: float = 0.0
fwhm_lp_bef: float = 0.0
imlib: Enum = 'vip-fft'
interpolation: Enum = 'lanczos4'
mask_rad: float = 0.0
mode: str = 'mean'
nproc: int = 1
verbose: bool = True
vip_hci.psfsub.rollsub.roll_sub(*all_args: List, **all_kwargs: dict)[source]

Perform roll-subtraction, followed by derotation and stacking of residual images.

Reference: [SCH98].

Parameters:
  • all_args (list, optional) – Positionnal arguments for the roll_sub algorithm. Full list of parameters below.

  • all_kwargs (dictionary, optional) – Mix of keyword arguments that can initialize a ROLL_SUB_Params and the optional rot_options dictionary (with keywords border_mode, mask_val, edge_blend, interp_zeros, ker; see docstrings of vip_hci.preproc.frame_rotate). Can also contain a ROLL_SUB_Params object/dictionary named algo_params.

  • cube (3d numpy ndarray, or tuple of two 3d numpy ndarray) – Input cube. Can also be a 4d array, with images obtained with the 1st and 2nd roll angle values are provided in the 1st and 2nd dimension of the 4d cube, respectively.

  • angle_list (1d numpy ndarray, or list/tuple of 2 elements) – Roll angles associated to each frame. Can also be a list/tuple of 2 elements. In the latter case, if input cube is 3D, it will assume that the first half of the frames are associated to the first roll angle, and the second half to the second roll angle.

  • mode ({'mean', 'median', 'individual'}, str optional) – If mode is set to ‘mean’ or ‘median’, only the mean/median frame of the image sequence obtained at the first roll angle is subtracted to the mean/media image from the sequence obtained with the second roll angle, and vice-versa. If mode is set to ‘individual’ a pair-wise subtraction of individual images obtained at each roll angle is performed, following the same order as in the input cube. To work, this mode requires the same number of images obtained with each roll angle.

  • imlib (Enum, see vip_hci.config.paramenum.Imlib) – See the documentation of vip_hci.preproc.frame_rotate.

  • interpolation (Enum, see vip_hci.config.paramenum.Interpolation) – See the documentation of the vip_hci.preproc.frame_rotate function.

  • collapse (Enum, see vip_hci.config.paramenum.Collapse) – Sets how temporal residual frames should be combined to produce an ADI image.

  • fwhm_lp_bef (float, optional) – FWHM of the Gaussian kernel used for low-pass filtering of the input images. Can be useful for mode=’individual’ and spatially subsampled input images. If set to 0 (default), no low-pass filtering is performed.

  • fwhm_lp_aft (float, optional) – FWHM of the Gaussian kernel used for low-pass filtering of the final image. Can be useful for spatially subsampled input images. If set to 0 (default), no low-pass filtering is performed.

  • nproc (None or int, optional) – Number of processes for parallel computing. If None the number of processes will be set to cpu_count()/2. By default the algorithm works in single-process mode.

  • full_output (bool, optional) – Whether to return the final median combined image only or with other intermediate arrays.

  • verbose (bool, optional) – If True prints to stdout intermediate info.

Returns:

  • cube_res (numpy ndarray, 3d) – [full_output=True] The cube of residuals.

  • cube_der (numpy ndarray, 3d) – [full_output=True] The derotated cube of residuals.

  • frame (numpy ndarray, 2d) – Median combination of the de-rotated cube.

vip_hci.psfsub.svd module

Module with functions for computing SVDs.

class vip_hci.psfsub.svd.SVDecomposer(data, mode='fullfr', inrad=None, outrad=None, svd_mode='lapack', scaling='temp-standard', scale_list=None, verbose=True)[source]

Bases: object

Class for SVD decomposition of 2d, 3d or 4d HCI arrays.

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

  • mode ({'fullfr', 'annular'}, optional) – Whether to use the whole frames or a single annulus.

  • inrad (None or int, optional) – [mode=’annular’] Inner radius.

  • outrad (None or int, optional) – [mode=’annular’] Outer radius.

  • svd_mode ({'lapack', 'arpack', 'eigen', 'randsvd', 'cupy', 'eigencupy',) –

    ‘randcupy’, ‘pytorch’, ‘eigenpytorch’, ‘randpytorch’}, str optional 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), proposed in [HAL09].

    • 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 (through Cupy). `

    • pytorch: uses the Pytorch library for GPU computation of the SVD.

    • eigenpytorch: offers the same method as with the eigen option but on GPU (through Pytorch).

    • randpytorch: is an adaptation of the randomized_svd algorithm, where all the linear algebra computations are done on a GPU (through Pytorch).

  • scaling ({None, "temp-mean", spat-mean", "temp-standard",) –

    “spat-standard”}, None or str optional Pixel-wise scaling mode using sklearn.preprocessing.scale function. If set to None, the input matrix is left untouched. Otherwise:

    • temp-mean: temporal px-wise mean is subtracted.

    • spat-mean: spatial mean is subtracted.

    • temp-standard: temporal mean centering plus scaling pixel values to unit variance (temporally).

    • spat-standard: spatial mean centering plus scaling pixel values to unit variance (spatially).

    DISCLAIMER: Using temp-mean or temp-standard scaling can improve the speckle subtraction for ASDI or (A)RDI reductions. Nonetheless, this involves a sort of c-ADI preprocessing, which (i) can be dangerous for datasets with low amount of rotation (strong self-subtraction), and (ii) should probably be referred to as ARDI (i.e. not RDI stricto sensu).

  • scale_list (numpy ndarray, optional) – Scaling factors to be used for re-scaling the spectral channels and aligning the speckles in the case of 4D HCI cubes.

  • verbose (bool, optional) – If True intermediate messages and timing are printed.

cevr_to_ncomp(cevr=0.9)[source]

Infer the number of principal components for a given CEVR.

Parameters:

cevr (float or tuple of floats, optional) – The desired CEVR.

Returns:

ncomp – The found number(s) of PCs.

Return type:

int or list of ints

generate_matrix()[source]

Generate a matrix from the input data. Pixel values in the matrix are scaled. Depending on mode, the matrix can come from an annulus instead of the whole frames.

get_cevr(ncomp_list=None, plot=True, plot_save=False, plot_dpi=90, plot_truncation=None)[source]

Calculate the cumulative explained variance ratio for the SVD of a cube/matrix (either full frames or a single annulus could be used).

Parameters:
  • ncomp_list (None, list or tuple, optional) – If provided the list is used to filter the vector of CEVR.

  • plot (bool, optional) – If True, the CEVR is plotted.

  • plot_save (bool, optional) – If True, the plot is saved as ./figure.pdf.

  • plot_dpi (int, optional) – The DPI of the figure.

  • plot_truncation (None or int, optional) – If provided, it created a second panel in the plot, focusing on the CEVR curve up to plot_truncation components.

Returns:

  • df_allks (pandas dataframe) – [ncomp_list is None] A table with the explained varaince ratio and the CEVR for all ncomps.

  • df_klist (pandas dataframe) – [ncomp_list is not None] A table with the ncomp_list, the explained varaince ratio and the CEVR.

run()[source]

Decompose the input data.

vip_hci.psfsub.utils_pca module

Module with helping functions for PCA.

vip_hci.psfsub.utils_pca.pca_annulus(cube, angs, ncomp, annulus_width, r_guess, cube_ref=None, svd_mode='lapack', scaling=None, collapse='median', weights=None, collapse_ifs='mean', **rot_options)[source]

PCA-ADI or PCA-RDI processed only for an annulus of the cube, with a given width and at a given radial distance to the frame center. It returns a processed frame with non-zero values only at the location of the annulus.

Parameters:
  • cube (3d or 4d numpy ndarray) – Input data cube to be processed by PCA.

  • angs (numpy ndarray or None) – The parallactic angles expressed as a numpy.array.

  • ncomp (int or list/1d numpy array of int) – The number of principal component.

  • annulus_width (float) – The annulus width in pixel on which the PCA is performed.

  • r_guess (float) – Radius of the annulus in pixels.

  • cube_ref (3d or 4d numpy ndarray, or list of 3d numpy ndarray, optional) – Reference library cube for Reference Star Differential Imaging.

  • svd_mode ({'lapack', 'randsvd', 'eigen', 'arpack'}, str optional) – Switch for different ways of computing the SVD and selected PCs.

  • scaling ({None, "temp-mean", spat-mean", "temp-standard",) –

    “spat-standard”}, None or str optional Pixel-wise scaling mode using sklearn.preprocessing.scale function. If set to None, the input matrix is left untouched. Otherwise:

    • temp-mean: temporal px-wise mean is subtracted.

    • spat-mean: spatial mean is subtracted.

    • temp-standard: temporal mean centering plus scaling pixel values to unit variance (temporally).

    • spat-standard: spatial mean centering plus scaling pixel values to unit variance (spatially).

    DISCLAIMER: Using temp-mean or temp-standard scaling can improve the speckle subtraction for ASDI or (A)RDI reductions. Nonetheless, this involves a sort of c-ADI preprocessing, which (i) can be dangerous for datasets with low amount of rotation (strong self-subtraction), and (ii) should probably be referred to as ARDI (i.e. not RDI stricto sensu).

  • collapse ({'median', 'mean', 'sum', 'wmean'}, str or None, optional) – Sets the way of collapsing the residual frames to produce a final image. If None then the cube of residuals is returned.

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

  • collapse_ifs ({'median', 'mean', 'sum', 'wmean', 'absmean'}, str or None, optional) – Sets the way of collapsing the spectral frames for producing a final image (in the case of a 4D input cube).

  • rot_options (dictionary, optional) – Dictionary with optional keyword values for “nproc”, “imlib”, “interpolation, “border_mode”, “mask_val”, “edge_blend”, “interp_zeros”, “ker” (see documentation of vip_hci.preproc.frame_rotate)

Returns:

pca_res – Depending on collapse and angs parameters, either a final collapsed frame (collapse not None) or the cube of residuals (derotated if angs is not None, non-derotated otherwises). Note: for 4d input cube, collapse must be non-None.

Return type:

2d or 3d ndarray

vip_hci.psfsub.utils_pca.pca_grid(cube, angle_list, fwhm=None, range_pcs=None, source_xy=None, cube_ref=None, mode='fullfr', annulus_width=20, svd_mode='lapack', scaling=None, mask_center_px=None, fmerit='mean', collapse='median', ifs_collapse_range='all', verbose=True, full_output=False, debug=False, plot=True, save_plot=None, start_time=None, scale_list=None, initial_4dshape=None, weights=None, exclude_negative_lobes=False, **rot_options)[source]

Compute a grid, depending on range_pcs, of residual PCA frames out of a 3d ADI cube (or a reference cube). If source_xy is provided, the number of principal components are optimized by measuring the S/N at this location on the frame (ADI, RDI). The metric used, set by fmerit, could be the given pixel’s S/N, the maximum S/N in a FWHM circular aperture centered on the given coordinates or the mean S/N in the same circular aperture.

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

  • angle_list (numpy ndarray, 1d) – Corresponding parallactic angle for each frame.

  • fwhm (None or float, optional) – Size of the FWHM in pixels, used for computing S/Ns when source_xy is passed.

  • range_pcs (None or tuple or list, optional) – The interval of PCs to be tried. If range_pcs is a tuple entered as (PC_INI, PC_MAX) a sequential grid will be evaluated between PC_INI and PC_MAX with step of 1. If range_pcs is a tuple entered as (PC_INI, PC_MAX, STEP] a grid will be evaluated between PC_INI and PC_MAX with the given STEP. If range_pcs is None, PC_INI=1, PC_MAX=n_frames-1 and STEP=1, which will result in longer running time. If range_pcs is a list, it should correspond to the number of PCs to be tested.

  • source_xy (None or tuple of floats) – X and Y coordinates of the pixel where the source is located and whose SNR is going to be maximized.

  • cube_ref (numpy ndarray, 3d, optional) – Reference library cube. For Reference Star Differential Imaging.

  • mode ({'fullfr', 'annular'}, optional) – Mode for PCA processing (full-frame or just in an annulus). There is a catch: the optimal number of PCs in full-frame may not coincide with the one in annular mode. This is due to the fact that the annulus matrix is smaller (less noisy, probably not containing the central star) and also its intrinsic rank (smaller that in the full frame case).

  • annulus_width (float, optional) – Width in pixels of the annulus in the case of the “annular” mode.

  • svd_mode ({'lapack', 'arpack', 'eigen', 'randsvd', 'cupy', 'eigencupy',) –

    ‘randcupy’, ‘pytorch’, ‘eigenpytorch’, ‘randpytorch’}, str optional 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), proposed in [HAL09].

    • 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 (through Cupy). `

    • pytorch: uses the Pytorch library for GPU computation of the SVD.

    • eigenpytorch: offers the same method as with the eigen option but on GPU (through Pytorch).

    • randpytorch: is an adaptation of the randomized_svd algorithm, where all the linear algebra computations are done on a GPU (through Pytorch).

  • scaling ({None, "temp-mean", spat-mean", "temp-standard",) –

    “spat-standard”}, None or str optional Pixel-wise scaling mode using sklearn.preprocessing.scale function. If set to None, the input matrix is left untouched. Otherwise:

    • temp-mean: temporal px-wise mean is subtracted.

    • spat-mean: spatial mean is subtracted.

    • temp-standard: temporal mean centering plus scaling pixel values to unit variance (temporally).

    • spat-standard: spatial mean centering plus scaling pixel values to unit variance (spatially).

    DISCLAIMER: Using temp-mean or temp-standard scaling can improve the speckle subtraction for ASDI or (A)RDI reductions. Nonetheless, this involves a sort of c-ADI preprocessing, which (i) can be dangerous for datasets with low amount of rotation (strong self-subtraction), and (ii) should probably be referred to as ARDI (i.e. not RDI stricto sensu).

  • mask_center_px (None or int, optional) – If None, no masking is done. If an integer > 1 then this value is the radius of the circular mask.

  • fmerit ({'px', 'max', 'mean'}) – The function of merit to be maximized. ‘px’ is source_xy pixel’s SNR, ‘max’ the maximum SNR in a FWHM circular aperture centered on source_xy and ‘mean’ is the mean SNR in the same circular aperture.

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

  • ifs_collapse_range (str 'all' or tuple of 2 int) – If a tuple, it should contain the first and last channels where the mSDI residual channels will be collapsed (by default collapses all channels).

  • verbose (bool, optional) – If True prints intermediate info and timing.

  • full_output (bool, optional) – If True it returns the optimal number of PCs, the final PCA frame for the optimal PCs and a cube with all the final frames for each number of PC that was tried.

  • debug (bool, bool optional) – Whether to print debug information or not.

  • plot (bool, optional) – Whether to plot the SNR and flux as functions of PCs and final PCA frame or not.

  • save_plot (string) – If provided, the pc optimization plot will be saved to that path.

  • start_time (None or datetime.datetime, optional) – Used when embedding this function in the main pca function. The object datetime.datetime is the global starting time. If None, it initiates its own counter.

  • scale_list (None or numpy ndarray, optional) – Scaling factors in case of IFS data (ADI+mSDI cube). They will be used to descale the spectral cubes and obtain the right residual frames, assuming cube is a 4d ADI+mSDI cube turned into 3d.

  • initial_4dshape (None or tuple, optional) – Shape of the initial ADI+mSDI cube.

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

  • exclude_negative_lobes (bool, opt) – Whether to include the adjacent aperture lobes to the tested location or not. Can be set to True if the image shows significant neg lobes.

  • rot_options (dictionary, optional) – Dictionary with optional keyword values for “nproc”, “imlib”, “interpolation”, “border_mode”, “mask_val”, “edge_blend”, “interp_zeros”, “ker” (see documentation of vip_hci.preproc.frame_rotate)

Returns:

  • cubeout (numpy ndarray) – 3D array with the residuals frames.

  • pclist (list) – [full_output=True, source_xy is None] The list of PCs at which the residual frames are computed.

  • finalfr (numpy ndarray) – [source_xy is not None] Residual frame with the highest S/N for source_xy.

  • df (pandas Dataframe) – [source_xy is not None] Dataframe with the pcs, measured fluxed and S/Ns for source_xy.

  • opt_npc (int) – [source_xy is not None] Optimal number of PCs for source_xy.

Module contents

Subpackage psfsub contains a large number of stellar PSF modelling + subtraction algorithms. The following methods have been implemented: - roll subtraction [SCH98] - median ADI/SDI [MAR06] / [SPA02] / [THA07] - frame differencing - a simplified version of LOCI [LAF07] - different flavours of PCA [AMA12] / [SOU12] working in full-frame, incremental and annular mode, including improvements, speed tricks, compatibility with ADI/RDI/SDI datasets and datasets too large to fit in memory. - full-frame and annular versions of NMF [LEE99] / [GOM17].