starred.deconvolution package
Submodules
starred.deconvolution.deconvolution module
- class starred.deconvolution.deconvolution.Deconv(image_size, number_of_sources, scale, upsampling_factor=2, epochs=1, psf=None, gaussian_fwhm=2, N_sersic=0, convolution_method='fft')[source]
Bases:
objectImage deconvolution class. The order of the params below matters for unpacking re-packing dictionaries.
- Parameters:
image_size (int) – input image size in pixels
number_of_sources (int) – amount of point sources
scale (float) – point source scaling factor
upsampling_factor (int) – the proportion by which the sampling rate increases with respect to the input image
epochs (int) – for a joint deconvolution, the number of epochs considered. For single deconvolution,
epochsis 1psf – Point Spread Function array from
starred.psf.psf, acting as deconvolution kernel heregaussian_fwhm (int) – the Gaussian’s FWHM in pixels
convolution_method (str) – method to use to calculate the convolution : ‘fft’, ‘scipy’ (recommended), or ‘lax’
- background_interpolate(kwargs, interpolation_order=0, inner_mask_size=3, outer_mask_size=7)[source]
Return an interpolated version of the background, interpolated below the point sources.
- Parameters:
kwargs – dictionary containing all keyword arguments
interpolation_order – order of the interpolation. For now, only 0 is available. It will take the mean of the pixel in the outer mask.
inner_mask_size – size of the region to interpolate
outer_mask_size – size of the outer mask used for the interpolation
- dump(path, kwargs, data=None, sigma_2=None, save_output_level=4, format='hdf5')[source]
Stores the model in a given file in pickle or hdf5 format (recommended).
- Parameters:
path – Filename of the output.
kwargs – Dictionary containing the fitted value of the model
data – (Nepoch x image_size x image_size) array containing the data
sigma_2 – (Nepoch x image_size x image_size) array containing the noise maps
save_output_level – Int. Level of output product to save: 1-just the parameters of the model, 2- add the input data, 3- add the output products (background, and pts source channel) 4- add the output products for every image.
- export(output_folder, kwargs_final, data, sigma_2, epoch=None, norm=1, format='fits', header=None)[source]
Saves all the output files.
- Parameters:
output_folder (str) – path to the output folder
kwargs_final – dictionary containing all keyword arguments
data – array containing the images
sigma_2 – array containing the noise maps
epoch – integer or array or list containing the indices of the epoch to export. None, for exporting all epochs.
norm – normalisation, to have the outputs in the correct unit
format (str) – output format. Choose between
npyorfitsheader – HDU header, as return by astropy.io.fits (optional). This is used to propagate the WCS coordinate to the outputfiles
- flux_at_epoch(kwargs, epoch=0)[source]
Return an array containing the flux of the point sources at epoch epoch.
- Parameters:
kwargs – dictionary containing all keyword arguments
epoch – index of the epoch
- getDeconvolved(kwargs, epoch=None)[source]
This is a utility function to retrieve the deconvolved model at one specific epoch. It will typically be called only once after the optimization process.
- Parameters:
kwargs – dictionary containing the parameters of the model
epoch (int) – number of epochs, if none, it will return the mean image
- Returns:
tuple of two 2D arrays: (deconvolved image, background)
- getPts_source(kwargs, epoch=None, source_ID=None)[source]
This is a utility function to retrieve an image with only the point source components at one specific epoch.
- Parameters:
kwargs – dictionary containing the parameters of the model
epoch (int) – number of epochs. If none, it will return the mean image
- Returns:
2D array.
- get_mask_pts_source(kwargs, mask_size=3, nan_mask=False)[source]
Return an array containing a mask of size ‘mask_size’ around the point sources. Output array is at the sub-sampled resolution.
- Parameters:
kwargs – dictionary containing all keyword arguments
mask_size – size of the mask around the point sources. In “small” pixel.
nan_mask – If true, return a mask where 0. have been replaced by nan
- mean_pts_amps(a)[source]
Return the mean amplitude per points source, over all epoch. Return a 1D vector of size self.M.
- Parameters:
a – amplitude vector, correspond to kwargs[‘kwargs_analytic’][‘a’]
- model(kwargs)[source]
Creates the 2D deconvolution model image.
- Parameters:
kwargs – dictionary containing the parameters of the model.
- Returns:
3D array - A stack of 2D images (one layer per epoch)
- modelstack(pars, h)[source]
Creates a model for all the epochs with shape (epoch, image_size_up, image_size_up) by vectorizing
self.one_epoch. Much faster for a large number of epochs (> 10), regime where theself.modelstack_forloopversion cannot be jit’ed anymore.- Parameters:
pars – stack of arrays of parameters for one epoch. Shape (epoch, N), where N is the number of parameters per epoch (see the docstring of
self.one_epoch)h – 1D array containing the pixelated background. Separated from pars because does not vary between epochs.
- Returns:
2D array of shape (epoch, image_size_up, image_size_up)
- modelstack_forloop(pars, h)[source]
Creates a model for all the epochs with shape (epoch, image_size_up, image_size_up) by looping over the epochs. This is slightly faster for a low number of epochs (~1 to 5) than using a vectorized version of
self.one_epoch.- Parameters:
pars – stack of arrays of parameters for one epoch. Shape (epoch, N), where N is the number of parameters per epoch (see the docstring of
self.one_epoch)h – 1D array containing the pixelated background. Separated from pars because does not vary between epochs.
- Returns:
2D array of shape (epoch, image_size_up, image_size_up)
- one_epoch(psf, params, h)[source]
Produces a 2D array corresponding to one epoch of the deconvolution.
- Parameters:
psf – 2D array representing the PSF at this epoch, usually upsampled
params –
1D array with the parameters in the following order:
the x positions of the point sources (M params)
the y positions of the point sources (M params)
the amplitudes of the point sources (M params)
the x translation of the whole epoch (1 param)
the y translation of the whole epoch (1 param)
h – a 1D array containing the unique background of the deconvolution,
self.image_size_up**2params.
- Returns:
2D array that represents the model for one epoch. This is, a collection of point sources plus
background, convolved with the PSF and potentially translated in the x and y directions
- sersic_background(R_sersic=[1.0], n_sersic=[1.0], e1=[0.0], e2=[0.0], center_x=[0.0], center_y=[0.0], amp=[0.0])[source]
Return a 1D vector containing all the pixels of the Sersic background.
- Parameters:
amp – surface brightness/amplitude value at the half light radius
R_sersic – semi-major axis half light radius, in big pixel
n_sersic – Sersic index
e1 – eccentricity parameter
e2 – eccentricity parameter
center_x – center in x-coordinate
center_y – center in y-coordinate
- setUpsamplingFactor(upsampling_factor)[source]
Sets the upsampling factor.
- Parameters:
upsampling_factor (int) – the proportion by which the sampling rate increases with respect to the input image
- shifted_gaussians(c_x, c_y, a, source_ID=None)[source]
Generates a 2D array with the point sources of the deconvolved image.
- Parameters:
c_x – 1D array containing the x positions of the point sources in unit of “big” pixel
c_y – 1D array containing the y positions of the point sources in unit of “big” pixel
a – 1D array containing the amplitude coefficients of each Gaussian representing a point source
- Returns:
2D array
- starred.deconvolution.deconvolution.load_Deconv_model(input_file, format='hdf5')[source]
Load Deconv model class from hdf5 or pickle file
- starred.deconvolution.deconvolution.setup_model(data, sigma_2, s, xs, ys, subsampling_factor, initial_a=None, astrometric_bound=None, dithering_bound=None, convolution_method='scipy', N_sersic=0, rotate=False)[source]
Utility setting up a deconvolution model. The returned dictionaries of parameters can later be adjusted by the user.
- Parameters:
data – 3D array containing the images, one per epoch. shape (epochs, im_size, im_size)
sigma_2 – 3D array containing the noisemaps, one per epoch. shape (epochs, im_size, im_size)
s – 3D array containing the narrow PSFs, one per epoch. shape (epochs, im_size_up, im_size_up) where im_size_up needs be a multiple of im_size.
xs – 1D array or list containing the x positions of the point sources. For M point sources, len(xs) is M. In units of big pixel.
ys – 1D array or list containing the y positions of the point sources. For M point sources, len(ys) is M. In units of big pixel.
initial_a – list containing the amplitudes of the point sources. For M point sources and N epochs, len(initial_a) is M*N. If none, amplitudes are applied scaled by the maax of each epoch.
subsampling_factor – integer, ratio of the size of the data pixels to that of the PSF pixels.
astrometric_bound – integer, maximum shift of the point sources, relative to the initial position. None: no upper nor lower limit
dithering_bound – integer, maximum shift from image to image. None: no upper nor lower limit
convolution_method – ‘scipy’, ‘fft’, or ‘lax’. To be passed to the Deconv class. Recommended : ‘scipy’.
N_sersic – Number of sersic profile to add to the background in addition to the pixelated grid (Default: 0)
rotate – bool, if true, takes rotation between epochs into account. If false, alpha gets fixed.
- Returns:
a tuple: (a starred.deconvolution.Deconv instance, a dictionary with the initial values of the model parameters, a dictionary with the upper boundaries of the model parameters, a dictionary with the lower boundaries of the model parameters, a dictionary with the fixed parameters of the model)
starred.deconvolution.loss module
- class starred.deconvolution.loss.Loss(data, deconv_class, param_class, sigma_2, regularization_terms='l1_starlet', regularization_strength_scales=1.0, regularization_strength_hf=1.0, regularization_strength_pts_source=0.0, regularization_strength_positivity=0.0, regularization_strength_positivity_ps=0.0, regularization_strength_flux_uniformity=0.0, regularize_full_model=False, W=None, prior=None)[source]
Bases:
objectClass that manages the (auto-differentiable) loss function, defined as: L = - log(likelihood) - log(regularization) - log(positivity) - log(prior)
Note that gradient, hessian, etc. are computed in the
InferenceBaseclass.- Parameters:
data – array containing the observations
deconv_class – deconvolution class from
starred.deconvolution.deconvolutionparam_class – parameters class from
starred.deconvolution.parameterssigma_2 – array containing the square of the noise maps
N (int) – number of observations stamps
regularization_terms (str) – Choose the type of regularization ‘l1_starlet’ or ‘l1_starlet_pts_rw’
regularization_strength_scales (float) – Lagrange parameter that weights intermediate scales in the transformed domain
regularization_strength_hf (float) – Lagrange parameter weighting the highest frequency scale
regularization_strength_positivity (float) – Lagrange parameter weighting the positivity of the background. 0 means no positivity constrain.
regularization_strength_positivity_ps – Lagrange parameter weighting the positivity of the point sources. 0 means no positivity constrain.
regularization_strength_pts_source (float) – Lagrange parameter regularising the point source channel.
regularization_strength_flux_uniformity – Lagrange parameter regularising scatter in the fluxes of each point source
W (jax.numpy.array) – weight matrix. Shape (n_scale, n_pix*subsampling_factor, n_pix*subsampling_factor)
regularize_full_model (bool) – option to regularise just the background (False) or background + point source channel (True)
prior (Prior class) – Prior class containing Gaussian prior on the free parameters
- property data
Returns the observations array.
- reduced_chi2(kwargs)[source]
Return the reduced chi2, given some model parameters
- Parameters:
kwargs – dictionary containing all keyword arguments
- property sigma_2
Returns the noise map array.
- class starred.deconvolution.loss.Prior(prior_analytic=None, prior_background=None, prior_sersic=None)[source]
Bases:
object- Parameters:
prior_analytic – list of [param_name, mean, 1-sigma priors]
prior_background – list of [param_name, mean, 1-sigma priors]
prior_sersic – list of [param_name, mean, 1-sigma priors]
starred.deconvolution.parameters module
- class starred.deconvolution.parameters.ParametersDeconv(kwargs_init, kwargs_fixed, kwargs_up=None, kwargs_down=None)[source]
Bases:
ParametersDeconvolution parameters class.
- Parameters:
kwargs_init – dictionary with information on the initial values of the parameters
kwargs_fixed – dictionary containing the fixed parameters
kwargs_up – dictionary with information on the upper bounds of the parameters
kwargs_down – dictionary with information on the lower bounds of the parameters
- get_param_names_for_model(kwargs_key)[source]
Returns the names of the parameters according to the key provided.
- kwargs2args(kwargs)[source]
Obtain an array of positional arguments from a dictionary of keyword arguments.
- param_names_analytic = ['a', 'c_x', 'c_y', 'dx', 'dy', 'alpha']
- param_names_background = ['mean', 'h']
- param_names_sersic = ['amp', 'R_sersic', 'n_sersic', 'center_x', 'center_y', 'e1', 'e2']
Module contents
This subpackage contains the deconvolution class/es