vip_hci.var package
Submodules
vip_hci.var.coords module
Module with functions related to image coordinates and coordinate conversions.
- vip_hci.var.coords.QU_to_QUphi(Q, U, delta_x=0, delta_y=0, scale_r2=False, north_convention=False)[source]
Returns Qphi and Uphi images, from input Q and U images.
- Parameters:
Q (numpy ndarray) – 2d numpy array containing the Q component of polarisation.
U (numpy ndarray) – 2d numpy array containing the U component of polarisation. Should have the same dimensions as Q.
delta_x (float, opt) – If the star is not at the center of the image, delta_x and delta_y indicate by how much it is offset along the x and y dimensions, resp.
delta_y (float, opt) – If the star is not at the center of the image, delta_x and delta_y indicate by how much it is offset along the x and y dimensions, resp.
scale_r2 (bool, opt) – Whether to scale by r^2 during conversion.
north_convention (bool, opt) – Whether to use angles measured from North up/East left (True), or measured from the positive x axis (False).
- Returns:
Qphi, Uphi – Qphi and Uphi images
- Return type:
numpy ndarrays
- vip_hci.var.coords.cart_to_pol(x, y, cx=0, cy=0, astro_convention=False)[source]
Returns polar coordinates for input cartesian coordinates
- Parameters:
x (float or numpy ndarray) – x coordinates with respect to the center
y (float or numpy ndarray) – y coordinates with respect to the center
cx (float or numpy ndarray) – x, y coordinates of the center of the image to be considered for conversion to cartesian coordinates.
cy (float or numpy ndarray) – x, y coordinates of the center of the image to be considered for conversion to cartesian coordinates.
astro_convention (bool) – Whether to use angles measured from North up/East left (True), or measured from the positive x axis (False).
- Returns:
r, theta – radii and polar angles corresponding to the input x and y.
- Return type:
floats or numpy ndarrays
- vip_hci.var.coords.dist(yc, xc, y1, x1)[source]
Return the Euclidean distance between two points, or between an array of positions and a point.
- vip_hci.var.coords.dist_matrix(n, cx=None, cy=None)[source]
Create matrix with euclidian distances from a reference point (cx, cy).
- Parameters:
n (int or 2D nd array) – output image shape is (n, n)
cx (float) – reference point. Defaults to the center.
cy (float) – reference point. Defaults to the center.
- Returns:
im
- Return type:
ndarray with shape (n, n)
Note
This is a replacement for ANDROMEDA’s DISTC.
- vip_hci.var.coords.frame_center(array, verbose=False)[source]
Return the coordinates y,x of the frame(s) center. If odd: dim/2-0.5 If even: dim/2
- Parameters:
array (2d/3d/4d numpy ndarray) – Frame or cube.
verbose (bool optional) – If True the center coordinates are printed out.
- Returns:
cy, cx – Coordinates of the center.
- Return type:
int
- vip_hci.var.coords.pol_to_cart(r, theta, r_err=0, theta_err=0, cx=0, cy=0, astro_convention=False)[source]
Returns cartesian coordinates for input polar coordinates, with error propagation.
- Parameters:
r (float or numpy ndarray) – radii and position angles to be converted to cartesian coords x and y.
theta (float or numpy ndarray) – radii and position angles to be converted to cartesian coords x and y.
r_err (float, optional) – Error on radial separation. Default is 0
theta_err (float, optional) – Error on position angle, in degrees. Default is 0
cx (float or numpy ndarray) – x, y coordinates of the center to be considered for conversion to cartesian coordinates.
cy (float or numpy ndarray) – x, y coordinates of the center to be considered for conversion to cartesian coordinates.
astro_convention (bool) – Whether to use angles measured from North up/East left (True), or measured from the positive x axis (False). If True, the x axis is reversed to match positive axis pointing East (left).
- Returns:
x, y (floats or numpy ndarrays) – x, y positions corresponding to input radii and position angles.
dx, dy (floats or numpy arrays) – dx, dy uncertainties on positions propagated from input uncertainties on r and theta.
- vip_hci.var.coords.pol_to_eq(r, t, rError=0, tError=0, astro_convention=False, plot=False)[source]
Converts a position (r,t) given in polar coordinates into \(\Delta\) RA and \(\Delta\) DEC (equatorial coordinates), with error propagation. Note: regardless of the assumption on input angle t (see description for astro_convention), the output RA is counted positive towards left.
- Parameters:
r (float) – The radial coordinate.
t (float) – The angular coordinate in degrees
rError (float, optional) – The error bar related to r.
tError (float, optional) – The error bar related to t, in deg.
astro_convention (bool, optional) – Whether the input angle t is assumed to be measured from North up, East left (True), or measured from the positive x axis (False).
plot (boolean, optional) – If True, a figure illustrating the error ellipse is displayed.
- Returns:
out – ((RA, RA error), (DEC, DEC error))
- Return type:
tuple
vip_hci.var.filters module
Module with frame/cube filtering functionalities.
- vip_hci.var.filters.cube_filter_highpass(array, mode='laplacian', verbose=True, **kwargs)[source]
Apply
frame_filter_highpass
to the frames of a 3d or 4d cube.- Parameters:
array (numpy ndarray) – Input cube, 3d or 4d.
mode (str, optional) –
mode
parameter to theframe_filter_highpass
function. Defaults to a Laplacian high-pass filter.verbose (bool, optional) – If
True
timing and progress bar are shown.**kwargs (dict) – Passed through to the
frame_filter_highpass
function.
- Returns:
filtered – High-pass filtered cube.
- Return type:
numpy ndarray
- vip_hci.var.filters.cube_filter_iuwt(cube, coeff=5, rel_coeff=1, full_output=False)[source]
Isotropic Undecimated Wavelet Transform filtering, as implemented in [KEN15] and detailed in [DAB15].
- Parameters:
cube (numpy ndarray) – Input cube.
coeff (int, optional) – Number of wavelet scales to be used in the decomposition.
rel_coeff (int, optional) – Number of relevant coefficients. In other words how many wavelet scales will represent in a better way our data. One or two scales are enough for filtering our images.
full_output (bool, optional) – If True, an additional cube with the multiscale decomposition of each frame will be returned.
- Returns:
cubeout (numpy ndarray) – Output cube with the filtered frames.
If full_output is True the filtered cube is returned together with the a
4d cube containing the multiscale decomposition of each frame.
- vip_hci.var.filters.cube_filter_lowpass(array, mode='gauss', median_size=5, fwhm_size=5, conv_mode='conv', kernel_sz=None, verbose=True, psf=None, mask=None, iterate=True, **kwargs)[source]
Apply
frame_filter_lowpass
to the frames of a 3d or 4d cube.- Parameters:
array (numpy ndarray) – Input cube, 3d or 4d.
mode (str, optional) – See the documentation of the
frame_filter_lowpass
function.median_size (int, optional) – See the documentation of the
frame_filter_lowpass
function.fwhm_size (float, optional) – See the documentation of the
frame_filter_lowpass
function.conv_mode (str, optional) – See the documentation of the
frame_filter_lowpass
function.kernel_sz (int, optional) – See the documentation of the
frame_filter_lowpass
function.verbose (bool, optional) – If True timing and progress bar are shown.
psf (numpy ndarray, optional) – Input normalised and centered psf, 2d frame. Should be provided if mode is set to ‘psf’.
mask (numpy ndarray, optional) – Binary mask indicating where the low-pass filtered image should be interpolated with astropy.convolution. This otion can be useful if the low-pass filtered image is aimed to capture low-spatial frequency sky signal, while avoiding a stellar halo (set to one in the binary mask). Note: only works with Gaussian kernel or PSF convolution.
iterate (bool, opt) – If the first convolution leaves nans, whether to continue replacing nans by interpolation until they are all replaced.
**kwargs (dict) – Passed through to the astropy.convolution.convolve or convolve_fft function.
- Returns:
filtered – Low-pass filtered cube.
- Return type:
numpy ndarray
- vip_hci.var.filters.frame_deconvolution(array, psf, n_it=30)[source]
Iterative image deconvolution following the scikit-image implementation of the Richardson-Lucy algorithm, described in [RIC72] and [LUC74].
Considering an image that has been convolved by the point spread function of an instrument, the algorithm will sharpen the blurred image through a user-defined number of iterations, which changes the regularisation.
- Parameters:
array (numpy ndarray) – Input image, 2d frame.
psf (numpy ndarray) – Input psf, 2d frame.
n_it (int, optional) – Number of iterations.
- Returns:
deconv – Deconvolved image.
- Return type:
numpy ndarray
- vip_hci.var.filters.frame_filter_highpass(array, mode, median_size=5, kernel_size=5, fwhm_size=5, btw_cutoff=0.2, btw_order=2, hann_cutoff=5, psf=None, conv_mode='conv', mask=None)[source]
High-pass filtering of input frame depending on parameter
mode
.The filtered image properties will depend on the
mode
and the relevant parameters.- Parameters:
array (numpy ndarray) – Input array, 2d frame.
mode (str) –
Type of High-pass filtering.
laplacian
applies a Laplacian filter with kernel size defined by
kernel_size
using the Opencv library.laplacian-conv
applies a Laplacian high-pass filter by defining a kernel (with
kernel_size
) and using theconvolve_fft
Astropy function.median-subt
subtracts a median low-pass filtered version of the image.
gauss-subt
subtracts a Gaussian low-pass filtered version of the image.
fourier-butter
applies a high-pass 2D Butterworth filter in Fourier domain.
hann
uses a Hann window.
median_size (int, optional) – Size of the median box for the
median-subt
filter.kernel_size (int, optional) – Size of the Laplacian kernel used in
laplacian
mode. It must be an positive odd integer value.fwhm_size (float, optional) – Size of the Gaussian kernel used in
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.psf (numpy ndarray, optional) – Input normalised and centered psf, 2d frame. Should be provided if mode is set to ‘psf’.
conv_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.
mask (numpy ndarray, optional) – Binary mask indicating where the low-pass filtered image should be interpolated with astropy.convolution. This otion can be useful if the low-pass filtered image is aimed to capture low-spatial frequency sky signal, while avoiding a stellar halo (set to one in the binary mask). Note: only works with Gaussian kernel or PSF convolution.
- Returns:
filtered – High-pass filtered image.
- Return type:
numpy ndarray
- vip_hci.var.filters.frame_filter_lowpass(array, mode='gauss', median_size=5, fwhm_size=5, conv_mode='convfft', kernel_sz=None, psf=None, mask=None, iterate=True, half_res_y=False, **kwargs)[source]
Low-pass filtering of input frame depending on parameter
mode
.- Parameters:
array (numpy ndarray) – Input array, 2d frame.
mode ({'median', 'gauss', 'psf'}, str optional) – Type of low-pass filtering.
median_size (int, optional) – Size of the median box for filtering the low-pass median filter.
fwhm_size (float or tuple of 2 floats, optional) – Size of the Gaussian kernel for the low-pass Gaussian filter. If a tuple is provided, it should correspond to y and x kernel sizes, respectively.
conv_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.
kernel_sz (int or None, optional) – Size of the kernel in pixels for 2D Gaussian and Moffat convolutions. If None, astropy.convolution will automatically consider 8*radius kernel sizes.
psf (numpy ndarray, optional) – Input normalised and centered psf, 2d frame. Should be provided if mode is set to ‘psf’.
mask (numpy ndarray, optional) – Binary mask indicating where the low-pass filtered image should be interpolated with astropy.convolution. This option can be useful if the low-pass filtered image is used to capture low-spatial frequency sky signal, while avoiding a stellar halo (set to one in the binary mask). Note: only works with Gaussian kernel or PSF convolution.
iterate (bool, opt) – If the first convolution leaves nans, whether to continue replacing nans by interpolation until they are all replaced.
half_res_y (bool, {True,False}, optional) – Whether the input data has only half the angular resolution vertically compared to horizontally (e.g. the case for some IFUs); in other words there are always 2 rows of pixels with exactly the same values. If so, the kernel will also be squashed vertically by a factor 2. Only used if mode is ‘gauss’
**kwargs (dict) – Passed through to the astropy.convolution.convolve or convolve_fft function.
- Returns:
filtered – Low-pass filtered image.
- Return type:
numpy ndarray
vip_hci.var.fit_2d module
2d fitting and creation of synthetic PSFs.
- vip_hci.var.fit_2d.create_synth_psf(model='gauss', shape=(9, 9), amplitude=1, x_mean=None, y_mean=None, fwhm=4, theta=0, gamma=None, alpha=1.5, radius=None, msdi=False)[source]
Creates a synthetic 2d or 3d PSF with a 2d model: Airy disk, Gaussian or Moffat, depending on
model
.- Parameters:
model ({'gauss', 'moff', 'airy'}, str optional) – Model to be used to create the synthetic PSF.
shape (tuple of ints, optional) – Shape of the output 2d array.
amplitude (float, optional) – Value of the amplitude of the 2d distribution.
x_mean (float or None, optional) – Value of the centroid in X of the distributions: the mean of the Gaussian or the location of the maximum of the Moffat or Airy disk models. If None, the centroid is placed at the center of the array.
y_mean (float or None, optional) – Value of the centroid in Y of the distributions: the mean of the Gaussian or the location of the maximum of the Moffat or Airy disk models. If None, the centroid is placed at the center of the array.
fwhm (float, tuple of floats, list or np.ndarray, optional) – FWHM of the model in pixels. For the Gaussian case, it controls the standard deviation of the Gaussian. If a tuple is given, then the Gaussian will be elongated (fwhm in x, fwhm in y). For the Moffat, it is related to the gamma and alpha parameters. For the Airy disk, it is related to the radius (of the first zero) parameter. If
msdi
is True thenfwhm
must be a list of 1d np.ndarray (for example for SPHERE/IFS this sounds like a reasonable FWHM: np.linspace(4.5,6.7,39)).theta (float, optional) – Rotation angle in degrees of the Gaussian.
gamma (float or None, optional) – Gamma parameter of core width of the Moffat model. If None, then it is calculated to correspond to the given
fwhm
.alpha (float, optional) – Power index of the Moffat model.
radius (float or None, optional) – The radius of the Airy disk (radius of the first zero). If None, then it is calculated to correspond to the given
fwhm
.msdi (bool, optional) – Creates a 3d PSF, for emulating an IFS PSF.
- Returns:
im – 2d array with given
shape
and containing the synthetic PSF.- Return type:
numpy ndarray
Note
- More information can be found at the following links:
http://docs.astropy.org/en/stable/api/astropy.modeling.functional_models.Gaussian2D.html http://docs.astropy.org/en/stable/api/astropy.modeling.functional_models.Moffat2D.html http://docs.astropy.org/en/stable/api/astropy.modeling.functional_models.AiryDisk2D.html
https://www.gnu.org/software/gnuastro/manual/html_node/PSF.html web.ipac.caltech.edu/staff/fmasci/home/astro_refs/PSFtheory.pdf web.ipac.caltech.edu/staff/fmasci/home/astro_refs/PSFsAndSampling.pdf
- vip_hci.var.fit_2d.fit_2d2gaussian(array, crop=False, cent=None, cropsize=15, fwhm_neg=4, fwhm_pos=4, theta_neg=0, theta_pos=0, neg_amp=1, fix_neg=True, threshold=False, sigfactor=2, bpm=None, full_output=False, debug=True)[source]
Fitting a 2D superimposed double Gaussian (negative and positive) to the 2D distribution of the data (reproduce e.g. the effect of a coronagraph)
- Parameters:
array (numpy ndarray) – Input frame with a single PSF.
crop (bool, optional) – If True a square sub image will be cropped equal to cropsize.
cent (tuple of float, optional) – X,Y position of the source in the array for extracting the subimage. If None the center of the frame is used for cropping the subframe. If fix_neg is set to True, this will also be used as the fixed position of the negative gaussian.
cropsize (int, optional) – Size of the subimage.
fwhm_neg (float or tuple of floats, optional) – Initial values for the FWHM of the fitted negative and positive Gaussians, in px. If a tuple, should be the FWHM value along x and y.
fwhm_pos (float or tuple of floats, optional) – Initial values for the FWHM of the fitted negative and positive Gaussians, in px. If a tuple, should be the FWHM value along x and y.
theta_neg (float, optional) – Angle of inclination of the 2d Gaussian counting from the positive X axis (only matters if a tuple was provided for fwhm_neg or fwhm_pos).
theta_pos (float, optional) – Angle of inclination of the 2d Gaussian counting from the positive X axis (only matters if a tuple was provided for fwhm_neg or fwhm_pos).
neg_amp (float, optional) – First guess on the amplitude of the negative gaussian, relative to the amplitude of the positive gaussian (i.e. 1 means the negative gaussian has the same amplitude as the positive gaussian)
fix_neg (bool, optional) – Whether to fix the position and FWHM of the negative gaussian for a fit with less free parameters. In that case, the center of the negative gaussian is assumed to be cent
threshold (bool, optional) – If True the background pixels (estimated using sigma clipped statistics) will be replaced by small random Gaussian noise.
sigfactor (int, optional) – The background pixels will be thresholded before fitting a 2d Gaussian to the data using sigma clipped statistics. All values smaller than (MEDIAN + sigfactor*STDDEV) will be replaced by small random Gaussian noise.
bpm (2D numpy ndarray, optional) – Mask of bad pixels to not consider for the fit.
full_output (bool, optional) – If False it returns just the centroid, if True also returns the FWHM in X and Y (in pixels), the amplitude and the rotation angle, and the uncertainties on each parameter.
debug (bool, optional) – If True, the function prints out parameters of the fit and plots the data, model and residuals.
- Returns:
mean_y (float) – Source centroid y position on input array from fitting.
mean_x (float) – Source centroid x position on input array from fitting.
If
full_output
is True it returns a Pandas dataframe containing thefollowing columns
- for the positive gaussian
’amplitude’ (Float value. Amplitude of the Gaussian.)
’centroid_x’ (Float value. X coordinate of the centroid.)
’centroid_y’ (Float value. Y coordinate of the centroid.)
’fwhm_x’ (Float value. FWHM in X [px].)
’fwhm_y’ (Float value. FWHM in Y [px].)
’theta’ (Float value. Rotation angle of x axis)
- for the negative gaussian
’amplitude_neg’ (Float value. Amplitude of the Gaussian.)
’centroid_x_neg’ (Float value. X coordinate of the centroid.)
’centroid_y_neg’ (Float value. Y coordinate of the centroid.)
’fwhm_x_neg’ (Float value. FWHM in X [px].)
’fwhm_y_neg’ (Float value. FWHM in Y [px].)
’theta_neg’ (Float value. Rotation angle of x axis)
- vip_hci.var.fit_2d.fit_2dairydisk(array, crop=False, cent=None, cropsize=15, fwhm=4, threshold=False, sigfactor=6, bpm=None, full_output=True, debug=True)[source]
Fitting a 2D Airy to the 2D distribution of the data.
- Parameters:
array (numpy ndarray) – Input frame with a single PSF.
crop (bool, optional) – If True a square sub image will be cropped equal to cropsize.
cent (tuple of int, optional) – X,Y integer position of source in the array for extracting the subimage. If None the center of the frame is used for cropping the subframe (the PSF is assumed to be ~ at the center of the frame).
cropsize (int, optional) – Size of the subimage.
fwhm (float, optional) – Initial values for the FWHM of the fitted 2d Airy disk, in px.
threshold (bool, optional) – If True the background pixels (estimated using sigma clipped statistics) will be replaced by small random Gaussian noise.
sigfactor (int, optional) – The background pixels will be thresholded before fitting a 2d Moffat to the data using sigma clipped statistics. All values smaller than (MEDIAN + sigfactor*STDDEV) will be replaced by small random Gaussian noise.
bpm (2D numpy ndarray, optional) – Mask of bad pixels to not consider for the fit.
full_output (bool, optional) – If False it returns just the centroid, if True also returns the FWHM in X and Y (in pixels), the amplitude and the rotation angle.
debug (bool, optional) – If True, the function prints out parameters of the fit and plots the data, model and residuals.
- Returns:
mean_y (float) – Source centroid y position on input array from fitting.
mean_x (float) – Source centroid x position on input array from fitting.
If
full_output
is True it returns a Pandas dataframe containing thefollowing columns
’amplitude’ (Float value. Moffat Amplitude.)
’centroid_x’ (Float value. X coordinate of the centroid.)
’centroid_y’ (Float value. Y coordinate of the centroid.)
’fwhm’ (Float value. FWHM [px].)
- vip_hci.var.fit_2d.fit_2dgaussian(array, crop=False, cent=None, cropsize=15, fwhmx=4, fwhmy=4, theta=0, threshold=False, sigfactor=6, bpm=None, full_output=True, debug=True)[source]
Fitting a 2D Gaussian to the 2D distribution of the data.
- Parameters:
array (numpy ndarray) – Input frame with a single PSF.
crop (bool, optional) – If True a square sub image will be cropped equal to cropsize.
cent (tuple of int, optional) – X,Y integer position of source in the array for extracting the subimage. If None the center of the frame is used for cropping the subframe (the PSF is assumed to be ~ at the center of the frame).
cropsize (int, optional) – Size of the subimage.
fwhmx (float, optional) – Initial values for the standard deviation of the fitted Gaussian, in px.
fwhmy (float, optional) – Initial values for the standard deviation of the fitted Gaussian, in px.
theta (float, optional) – Angle of inclination of the 2d Gaussian counting from the positive X axis.
threshold (bool, optional) – If True the background pixels (estimated using sigma clipped statistics) will be replaced by small random Gaussian noise.
sigfactor (int, optional) – The background pixels will be thresholded before fitting a 2d Gaussian to the data using sigma clipped statistics. All values smaller than (MEDIAN + sigfactor*STDDEV) will be replaced by small random Gaussian noise.
bpm (2D numpy ndarray, optional) – Mask of bad pixels to not consider for the fit.
full_output (bool, optional) – If False it returns just the centroid, if True also returns the FWHM in X and Y (in pixels), the amplitude and the rotation angle, and the uncertainties on each parameter.
debug (bool, optional) – If True, the function prints out parameters of the fit and plots the data, model and residuals.
- Returns:
mean_y (float) – Source centroid y position on input array from fitting.
mean_x (float) – Source centroid x position on input array from fitting.
If
full_output
is True it returns a Pandas dataframe containing thefollowing columns – ‘centroid_y’: Y coordinate of the centroid. ‘centroid_x’: X coordinate of the centroid. ‘fwhm_y’: Float value. FWHM in X [px]. ‘fwhm_x’: Float value. FWHM in Y [px]. ‘amplitude’: Amplitude of the Gaussian. ‘theta’: Float value. Rotation angle. # and fit uncertainties on the above values: ‘centroid_y_err’ ‘centroid_x_err’ ‘fwhm_y_err’ ‘fwhm_x_err’ ‘amplitude_err’ ‘theta_err’
- vip_hci.var.fit_2d.fit_2dmoffat(array, crop=False, cent=None, cropsize=15, fwhm=4, threshold=False, sigfactor=6, bpm=None, full_output=True, debug=True)[source]
Fitting a 2D Moffat to the 2D distribution of the data.
- Parameters:
array (numpy ndarray) – Input frame with a single PSF.
crop (bool, optional) – If True a square sub image will be cropped equal to cropsize.
cent (tuple of int, optional) – X,Y integer position of source in the array for extracting the subimage. If None the center of the frame is used for cropping the subframe (the PSF is assumed to be ~ at the center of the frame).
cropsize (int, optional) – Size of the subimage.
fwhm (float, optional) – Initial values for the FWHM of the fitted 2d Moffat, in px.
threshold (bool, optional) – If True the background pixels (estimated using sigma clipped statistics) will be replaced by small random Gaussian noise.
sigfactor (int, optional) – The background pixels will be thresholded before fitting a 2d Moffat to the data using sigma clipped statistics. All values smaller than (MEDIAN + sigfactor*STDDEV) will be replaced by small random Gaussian noise.
bpm (2D numpy ndarray, optional) – Mask of bad pixels to not consider for the fit.
full_output (bool, optional) – If False it returns just the centroid, if True also returns the FWHM in X and Y (in pixels), the amplitude and the rotation angle.
debug (bool, optional) – If True, the function prints out parameters of the fit and plots the data, model and residuals.
- Returns:
mean_y (float) – Source centroid y position on input array from fitting.
mean_x (float) – Source centroid x position on input array from fitting.
If
full_output
is True it returns a Pandas dataframe containing thefollowing columns
’alpha’ (Float value. Alpha parameter.)
’amplitude’ (Float value. Moffat Amplitude.)
’centroid_x’ (Float value. X coordinate of the centroid.)
’centroid_y’ (Float value. Y coordinate of the centroid.)
’fwhm’ (Float value. FWHM [px].)
’gamma’ (Float value. Gamma parameter.)
vip_hci.var.iuwt module
Code with the Isotropic Undecimated Wavelet Transform, taken (Aug 24, 2015) from https://github.com/ratt-ru/PyMORESANE/ Credits to J. S. Kenyon
- vip_hci.var.iuwt.iuwt_decomposition(in1, scale_count, scale_adjust=0, mode='ser', core_count=2, store_smoothed=False)[source]
This function serves as a handler for the different implementations of the IUWT decomposition. It allows the different methods to be used almost interchangeably.
The code was taken from [KEN15] and is detailed in [DAB15].
INPUTS: in1 (no default): Array on which the decomposition is to be performed. scale_count (no default): Maximum scale to be considered. scale_adjust (default=0): Adjustment to scale value if first scales are of no interest. mode (default=’ser’): Implementation of the IUWT to be used - ‘ser’, ‘mp’. core_count (default=1): Additional option for multiprocessing - specifies core count. store_smoothed (default=False): Boolean specifier for whether the smoothed image is stored or not.
OUTPUTS: Returns the decomposition with the additional smoothed coefficients if specified.
- vip_hci.var.iuwt.iuwt_recomposition(in1, scale_adjust=0, mode='ser', core_count=1, store_on_gpu=False, smoothed_array=None)[source]
This function serves as a handler for the different implementations of the IUWT recomposition. It allows the different methods to be used almost interchangeably.
INPUTS: in1 (no default): Array on which the decomposition is to be performed. scale_adjust (no default): Number of omitted scales. mode (default=’ser’) Implementation of the IUWT to be used - ‘ser’, ‘mp’ or ‘gpu’. core_count (default=1) Additional option for multiprocessing - specifies core count. store_on_gpu (default=False): Boolean specifier for whether the decomposition is stored on the gpu or not.
OUTPUTS: Returns the recomposition.
- vip_hci.var.iuwt.mp_a_trous(C0, wavelet_filter, scale, core_count)[source]
This is a reimplementation of the a trous filter which makes use of multiprocessing. In particular, it divides the input array of dimensions NxN into M smaller arrays of dimensions (N/M)xN, where M is the number of cores which are to be used.
INPUTS: C0 (no default): The current array which is to be decomposed. wavelet_filter (no default): The filter-bank which is applied to the components of the transform. scale (no default): The scale at which decomposition is to be carried out. core_count (no default): The number of CPU cores over which the task should be divided.
OUTPUTS: shared_array The decomposed array.
- vip_hci.var.iuwt.mp_a_trous_kernel(C0, wavelet_filter, scale, slice_ind, slice_width, r_or_c='row')[source]
This is the convolution step of the a trous algorithm.
INPUTS: C0 (no default): The current array which is to be decomposed. wavelet_filter (no default): The filter-bank which is applied to the elements of the transform. scale (no default): The scale at which decomposition is to be carried out. slice_ind (no default): The index of the particular slice in which the decomposition is being performed. slice_width (no default): The number of elements of the shorter dimension of the array for decomposition. r_or_c (default = “row”): Indicates whether strips are rows or columns.
OUTPUTS: NONE - C0 is a mutable array which this function alters. No value is returned.
- vip_hci.var.iuwt.mp_iuwt_decomposition(in1, scale_count, scale_adjust, store_smoothed, core_count)[source]
This function calls the a trous algorithm code to decompose the input into its wavelet coefficients. This is the isotropic undecimated wavelet transform implemented for multiple CPU cores. NOTE: Python is not well suited to multiprocessing - this may not improve execution speed.
INPUTS: in1 (no default): Array on which the decomposition is to be performed. scale_count (no default): Maximum scale to be considered. scale_adjust (default=0): Adjustment to scale value if first scales are of no interest. store_smoothed (default=False):Boolean specifier for whether the smoothed image is stored or not. core_count (no default): Indicates the number of cores to be used.
OUTPUTS: detail_coeffs Array containing the detail coefficients. C0 (optional): Array containing the smoothest version of the input.
- vip_hci.var.iuwt.mp_iuwt_recomposition(in1, scale_adjust, core_count, smoothed_array)[source]
This function calls the a trous algorithm code to recompose the input into a single array. This is the implementation of the isotropic undecimated wavelet transform recomposition for multiple CPU cores.
INPUTS: in1 (no default): Array containing wavelet coefficients. scale_adjust (no default): Indicates the number of omitted array pages. core_count (no default): Indicates the number of cores to be used. smoothed_array (default=None): For a complete inverse transform, this must be the smoothest approximation.
OUTPUTS: recomposiiton Array containing the reconstructed image.
- vip_hci.var.iuwt.ser_a_trous(C0, filter, scale)[source]
The following is a serial implementation of the a trous algorithm. Accepts the following parameters:
INPUTS: filter (no default): The filter-bank which is applied to the components of the transform. C0 (no default): The current array on which filtering is to be performed. scale (no default): The scale for which the decomposition is being carried out.
OUTPUTS: C1 The result of applying the a trous algorithm to the input.
- vip_hci.var.iuwt.ser_iuwt_decomposition(in1, scale_count, scale_adjust, store_smoothed)[source]
This function calls the a trous algorithm code to decompose the input into its wavelet coefficients. This is the isotropic undecimated wavelet transform implemented for a single CPU core.
INPUTS: in1 (no default): Array on which the decomposition is to be performed. scale_count (no default): Maximum scale to be considered. scale_adjust (default=0): Adjustment to scale value if first scales are of no interest. store_smoothed (default=False):Boolean specifier for whether the smoothed image is stored or not.
OUTPUTS: detail_coeffs Array containing the detail coefficients. C0 (optional): Array containing the smoothest version of the input.
- vip_hci.var.iuwt.ser_iuwt_recomposition(in1, scale_adjust, smoothed_array)[source]
This function calls the a trous algorithm code to recompose the input into a single array. This is the implementation of the isotropic undecimated wavelet transform recomposition for a single CPU core.
INPUTS: in1 (no default): Array containing wavelet coefficients. scale_adjust (no default): Indicates the number of truncated array pages. smoothed_array (default=None): For a complete inverse transform, this must be the smoothest approximation.
OUTPUTS: recomposition Array containing the reconstructed image.
vip_hci.var.shapes module
Module with various functions to create shapes, annuli and segments.
- vip_hci.var.shapes.create_ringed_spider_mask(im_shape, ann_out, ann_in=0, sp_width=10, sp_angle=0, nlegs=6)[source]
Mask out information outside the annulus and inside the spiders (zeros).
- Parameters:
im_shape (tuple of int) – Tuple of length two with 2d array shape (Y,X).
ann_out (int) – Outer radius of the annulus.
ann_in (int, opt) – Inner radius of the annulus.
sp_width (int, opt) – Width of the spider arms (6 legs by default).
sp_angle (int, opt) – angle of the first spider arm (on the positive horizontal axis) in counter-clockwise sense.
nlegs (int, opt) – Number of legs of the spider.
- Returns:
mask – 2d array of zeros and ones.
- Return type:
numpy ndarray
- vip_hci.var.shapes.get_annular_wedge(data, inner_radius, width, wedge=(0, 360), mode='ind')[source]
Return indices or values in segments of a centered annulus.
The annulus is defined by
inner_radius <= annulus < inner_radius+width
.- Parameters:
data (2d numpy ndarray or tuple) – Input 2d array (image) ot tuple with its shape.
inner_radius (float) – The inner radius of the donut region.
width (float) – The size of the annulus.
wedge (tuple of 2 floats) – Initial and final azimuths [degrees] of the annular segment, counting from the positive x-axis counter-clockwise.
mode ({'ind', 'val', 'mask'}, optional) – Controls what is returned: indices of selected pixels, values of selected pixels, or a boolean mask.
- Returns:
indices (list of ndarrays) – [mode=’ind’] Coordinates of pixels for each annulus segment.
values (list of ndarrays) – [mode=’val’] Pixel values.
masked (list of ndarrays) – [mode=’mask’] Copy of
data
with masked out regions.
Note
Moving from
get_annulus
toget_annulus_segments
:
- vip_hci.var.shapes.get_annulus_segments(data, inner_radius, width, nsegm=1, theta_init=0, optim_scale_fact=1, mode='ind', out=False)[source]
Return indices or values in segments of a centered annulus.
The annulus is defined by
inner_radius <= annulus < inner_radius+width
.- Parameters:
data (2d numpy ndarray or tuple) – Input 2d array (image) ot tuple with its shape.
inner_radius (float) – The inner radius of the donut region.
width (float) – The size of the annulus.
nsegm (int) – Number of segments of annulus to be extracted.
theta_init (int) – Initial azimuth [degrees] of the first segment, counting from the positive x-axis counterclockwise.
optim_scale_fact (float) – To enlarge the width of the segments, which can then be used as optimization segments (e.g. in LOCI).
mode ({'ind', 'val', 'mask'}, optional) – Controls what is returned: indices of selected pixels, values of selected pixels, or a boolean mask.
out (bool; optional) – Return all indices or values outside the centered annulus.
- Returns:
indices (list of ndarrays) – [mode=’ind’] Coordinates of pixels for each annulus segment.
values (list of ndarrays) – [mode=’val’] Pixel values.
masked (list of ndarrays) – [mode=’mask’] Copy of
data
with masked out regions.
Note
Moving from
get_annulus
toget_annulus_segments
:
- vip_hci.var.shapes.get_circle(array, radius, cy=None, cx=None, mode='mask')[source]
Return a centered circular region from a 2d ndarray.
- Parameters:
array (numpy ndarray) – Input 2d array or image.
radius (int) – The radius of the circular region.
cy (int, optional) – Coordinates of the circle center. If one of them is
None
, the center ofarray
is used.cx (int, optional) – Coordinates of the circle center. If one of them is
None
, the center ofarray
is used.mode ({'mask', 'val'}, optional) – Controls what is returned: array with circular mask applied, or values of the pixels in the circular region.
- Returns:
masked (numpy ndarray) – [mode=”mask”] Input array with the circular mask applied.
values (numpy ndarray) – [mode=”val”] 1d array with the values of the pixels in the circular region.
Note
An alternative implementation would use
skimage.draw.disk
.disk
performs better on large arrays (e.g. 1000px, 10.000px), while the current implementation is faster for small arrays (e.g. 100px). See test_shapes.py for benchmark details.
- vip_hci.var.shapes.get_ell_annulus(data, a, b, PA, width, cy=None, cx=None, mode='ind')[source]
Return a centered elliptical annulus from a 2d ndarray.
All the rest pixels are set to zeros.
- Parameters:
data (numpy ndarray or tuple) – Input 2d array (image) or tuple with a shape.
a (float) – Semi-major axis.
b (float) – Semi-minor axis.
PA (deg, float) – The PA of the semi-major axis in degrees.
width (float) – The size of the annulus along the semi-major axis; it is proportionnally thinner along the semi-minor axis.
output_values ({False, True}, optional) – If True returns the values of the pixels in the annulus.
cy (int or None, optional) – Coordinates of the circle center. If
None
, the center is determined by theframe_center
function.cx (int or None, optional) – Coordinates of the circle center. If
None
, the center is determined by theframe_center
function.mode ({'ind', 'val', 'mask'}, optional) – Controls what is returned: indices of selected pixels, values of selected pixels, or a boolean mask.
- Returns:
indices (tuple(y, x)) – [mode=’ind’] Coordinates of the inner elliptical region.
values (1d ndarray) – [mode=’val’] Values of the pixels in the inner elliptical region.
masked (2d ndarray) – [mode=’mask’] Input image where the outer region is masked with
0
.
- vip_hci.var.shapes.get_ellipse(data, a, b, pa, cy=None, cx=None, mode='ind')[source]
Return a centered elliptical region from a 2d ndarray.
- Parameters:
data (numpy ndarray or tuple) – Input 2d array (image) or tuple with a shape.
a (float) – Semi-major axis.
b (float) – Semi-minor axis.
pa (deg, float) – The PA of the semi-major axis in degrees.
cy (int or None, optional) – Coordinates of the circle center. If
None
, the center is determined by theframe_center
function.cx (int or None, optional) – Coordinates of the circle center. If
None
, the center is determined by theframe_center
function.mode ({'ind', 'val', 'mask', 'bool'}, optional) – Controls what is returned: indices of selected pixels, values of selected pixels, or a boolean mask.
- Returns:
indices (tuple(y, x)) – [mode=’ind’] Coordinates of the inner elliptical region.
values (1d ndarray) – [mode=’val’] Values of the pixels in the inner elliptical region.
masked (2d ndarray) – [mode=’mask’] Input image where the outer region is masked with
0
.bool_mask (2d boolean ndarray) – [mode=’bool’] A boolean mask where
True
is the inner region.
- vip_hci.var.shapes.get_square(array, size, y, x, position=False, force=False, verbose=True)[source]
Return an square subframe from a 2d array or image.
- Parameters:
array (2d numpy ndarray) – Input frame.
size (int) – Size of the subframe.
y (int or float) – Y coordinate of the center of the subframe (obtained with the function
frame_center
).x (int or float) – X coordinate of the center of the subframe (obtained with the function
frame_center
).position (bool, optional) – If set to True return also the coordinates of the bottom-left vertex.
force (bool, optional) – Size and the size of the 2d array must be both even or odd. With
force
set to False, the requested size is flexible (i.e. +1 can be applied to requested crop size for its parity to match the input size). Ifforce
set to True, the requested crop size is enforced, even if parities do not match (warnings are raised!).verbose (bool optional) – If True, warning messages might be shown.
- Returns:
array_out (numpy ndarray) – Sub array.
y0, x0 (int) – [position=True] Coordinates of the bottom-left vertex.
- vip_hci.var.shapes.mask_circle(array, radius, fillwith=0, mode='in', cy=None, cx=None, output='masked_arr')[source]
Mask the pixels inside/outside of a centered circle with
fillwith
.Returns a modified copy of
array
.- Parameters:
array (2d/3d/4d numpy ndarray) – Input frame or cube.
radius (float) – Radius of the circular mask.
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 tofillwith
.cy (floats, opt) – XY coordinates of thenter of the mask. By default, it considers the center of the image.
cx (floats, opt) – XY coordinates of thenter of the mask. By default, it considers the center of the image.
output ({'masked_arr', 'bool_mask'}, optional) – Whether to return the masked frame or a bolean mask
- Returns:
array_masked – Masked frame or cube.
- Return type:
numpy ndarray
- vip_hci.var.shapes.mask_roi(array, source_xy, exc_radius=4, ann_width=4, inc_radius=8, mode='val', plot=False)[source]
Return a mask corresponding to the region of interest for a test point source as defined in [GEB20].
Given a frame and a location of interest with coordinates xy=(cx,cy), the mask consists of all pixels from three different regions with respect to xy:
(r1) Exclusion region: Pixels from the region of interest. These are excluded in the final mask.
(r2) Local region: Pixels around xy in a circular fashion.
(r3) Symmetric local region: Pixels around the (anti)symmetric xy with respect to the star location. It is also defined in a circular fashion with same radius as “local region.”
(r4) Annulus region: Pixels from the annulus where xy is located.
The goal of this mask is to disentangle the expected structure of the speckle pattern. [GEB20] comment that ‘r2 is chosen to capture any local effects around xy due to the instrument. r3 is chosen symmetrically opposite of xy because if there is a speckle at xy, there should also be an (anti-)speckle at r2. r4 is used because we know that the systematic noise significantly depends on the radial variable.’
- Parameters:
array (2d numpy ndarray) – Input 2d array (image) ot tuple with its shape.
source_xy (tuple of 2 int) – Pixel coordinates of the location of interest.
exc_radius (float) – Known size of the FWHM in pixels to be used. Default is 4.
ann_width (int) – Width (in pxs) of the annulus (region 3 on description).
inc_radius (int) – Radius (in pxs) of the circular regions r2 and r3 on description. Default is 8.
mode ({'bool', 'val', 'mask', 'ind'}, optional) – Controls what is returned: array with the mask applied, values of the pixels in the mask region or indices of selected pixels of the mask. Default is ‘mask’.
plot (bool, optional) – Whether to display a plot showing the defined mask.
- Returns:
bool_mask (numpy ndarray) – [mode=”bool”] Input array with the circular mask applied.
values (numpy ndarray) – [mode=”val”] 1d array with the values of the pixels in the circular region.
masked (numpy ndarray) – [mode=”mask”] Input array with the circular mask applied.
indices (tuple(y, x)) – [mode=’ind’] Coordinates of the mask.
- vip_hci.var.shapes.matrix_scaling(matrix, scaling)[source]
Scale a matrix using
sklearn.preprocessing.scale
function.- Parameters:
matrix (2d numpy ndarray) – Input 2d array.
scaling (None or string) –
Scaling method.
None
no scaling is performed on the input data before SVD
"temp-mean"
temporal px-wise mean subtraction
"spat-mean"
the spatial mean is subtracted
temp-standard"
temporal mean centering plus scaling to unit variance
"spat-standard"
spatial mean centering plus scaling to unit variance
- Returns:
matrix – 2d array with scaled values.
- Return type:
2d numpy ndarray
- vip_hci.var.shapes.prepare_matrix(array, scaling=None, mask_center_px=None, mode='fullfr', inner_radius=None, outer_radius=None, discard_mask_pix=False, verbose=True)[source]
Build the matrix for the SVD/PCA and other matrix decompositions.
Center the data and mask the frame’s central area if needed.
- Parameters:
array (3d numpy ndarray) – Input cube.
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
ortemp-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) – [mode=’fullfr’] Whether to mask the center of the frames or not.
mode ({'fullfr', 'annular'}, optional) – Whether to use the whole frames or a single annulus.
inner_radius (int or float, optional) – [mode=’annular’] Distance in pixels from the center of the frame to the inner radius of the annulus.
outer_radius (int or float, optional) – [mode=’annular’] Distance in pixels from the center of the frame to the outer radius of the annulus.
discard_mask_pix (bool, optional) – [mask_center_px=int ] Exclude masked pixels in the vetorized frames
verbose (bool, optional) – If True prints intermediate info.
- Returns:
matrix (2d numpy ndarray) – Out matrix whose rows are vectorized frames from the input cube.
ind (tuple) – [mode=’annular’] Indices of the annulus as
(yy, xx)
.
- vip_hci.var.shapes.reshape_matrix(array, y, x)[source]
Convert a matrix whose rows are vect. frames to a cube with reshaped frames.
- Parameters:
array (2d ndarray) – Input data of shape
(nframes, npixels)
. Every row (array[n]
) corresponds to one vectorized (“flattened”) 2d frame.y (int) – desired height and width of the frames.
y*x = npixels
x (int) – desired height and width of the frames.
y*x = npixels
- Returns:
cube – Cube of shape
(nframes, y, x)
.- Return type:
3d ndarray
Examples
In [1]: vect_frames = np.array([[1, 1, 1, 2, 2, 2], [1, 2, 3, 4, 5, 6]]) In [2]: cube = vip.var.reshape_matrix(vect_frames, 2, 3) In [3]: cube Out[3]: array([[[1, 1, 1], [2, 2, 2]], [[1, 2, 3], [4, 5, 6]]]) In [4]: cube.shape Out[4]: (2, 2, 3)
Module contents
Subpackage var
has helping functions such as:
image filtering,
shapes extraction (annulus, squares subimages, circular apertures),
2d fitting (Gaussian, Moffat),
wavelet transform definition,
Q/U to Qphi/Uphi image conversion (for polarimetric data).