mangadap.spectra.rowstackedspectra module

Base class for row-stacked spectra


License

Copyright © 2019, SDSS-IV/MaNGA Pipeline Group


class mangadap.spectra.rowstackedspectra.RowStackedSpectra(wave, flux, ivar=None, mask=None, bitmask=None, sres=None, xpos=None, ypos=None, area=None, log=True, prihdr=None, fluxhdr=None)[source]

Bases: object

Base container class for a set of 1D extracted spectra.

Row-stacked spectra have two axes: The first is just used to organize the spectra, whereas the second is the spectral axis. That is, one can loop through the spectra by, e.g.:

for s in spectra:
    plt.plot(s)
    ...

The spectral axis is expected to be the same for all spectra.

Currently no spectral or “spatial” (e.g. cross-talk) covariance is accounted for.

Warning

The class can be used to construct datacubes from row-stacked spectra. In these calculations, the aperture area is used to rescale the flux to ensure that it is conserved per unit area. While the class allows the aperture area to be fiber-specific, it does not yet allow the smoothing kernel used in the rectification to also be fiber-specific. Beware if you’re trying to construct a cube using multiple fiber formats!

Derived classes should, in particular, provide the read methods. The critical components of the derived classes are:

Todo

  • Fill in these “critical components”

  • How do I abstract how to provide the on-sky aperture?

Parameters:
  • wave (numpy.ndarray) – The wavelength vector of all spectra. Shape is \((N_\lambda,)\).

  • flux (numpy.ndarray) – Data array with the flux for each spectrum. Shape should be \((N_{\rm spec}, N_\lambda)\).

  • ivar (numpy.ndarray, optional) – The inverse variance of the flux array. Shape must be the same as flux. If None, errors are ignored.

  • mask (numpy.ndarray, optional) – Mask array for the flux measurements. Can be a boolean mask (False=unmasked,good; True=masked,bad) or an integer array associated with the provided mangadap.util.bitmask.BitMask object (see bitmask). Shape must be the same as flux.

  • bitmask (mangadap.util.bitmask.BitMask, optional) – Object used to select and toggle masked pixels from mask. If None, any provided mask must be a boolean array or can be successfully converted to one.

  • sres (numpy.ndarray, optional) – The spectral resolution, \(R = \lambda / \Delta\lambda\). Can be a single vector with shape \((N_\lambda,)\), if the spectral resolution is independent of spatial position, or a 2D set of spectral resolution vectors with a shape that matches flux. If None, spectral resolution is ignored.

  • xpos (numpy.ndarray, optional) – The wavelength-dependent coordinate along the right-ascension direction of the on-sky aperture for all spectra. Expected to be arcseconds from some fiducial coordinate. Coordinates should increase from W to E. Shape must be the same as flux.

  • ypos (numpy.ndarray, optional) – The wavelength-dependent coordinate along the declination direction of the on-sky aperture for all spectra. Expected to be arcseconds from some fiducial coordinate. Coordinates should increase from S to N. Shape must be the same as flux.

  • area (float, numpy.ndarray, optional) – The fiducial area subtended by each spectral aperture, assumed to be wavelength independent. The aperture is assumed to be circular. Can provide a single value for all fibers or an array with shape \((N_{\rm spec},)\) that provides fiber-specific values. If None, the area is set to unity, which has follow-on effects on the results provided by some of the class methods, as described where relevant.

  • log (bool, optional) – Flag that the datacube spectral pixels are binned logarithmically in wavelength.

  • prihdr (astropy.io.fits.Header, optional) – Primary header read from source fits file. If None, instantiated as an empty astropy.io.fits.Header.

  • fluxhdr (astropy.io.fits.Header, optional) – Header specifically for the flux extension of the source fits file. If None, set to be a copy of the primary header.

shape

Shape of the main spectral arrays.

Type:

tuple

nwave

Number of wavelength channels.

Type:

int

bitmask

Object used to select and toggle masked pixels from mask. Can be None.

Type:

mangadap.util.bitmask.BitMask

log

Flag that the spectral pixels are binned logarithmically in wavelength.

Type:

bool

wave

Wavelength vector applicable to all spectra.

Type:

numpy.ndarray

flux

Flux array.

Type:

numpy.ndarray

ivar

Flux inverse variance. Can be None.

Type:

numpy.ndarray

mask

Spectral mask. Can be None.

Type:

numpy.ndarray

sres

Spectral-resolution array. Can be None.

Type:

numpy.ndarray

xpos

The wavelength-dependent coordinate along the right-ascension direction of the on-sky aperture for all spectra. Expected to be arcseconds from some fiducial coordinate.

Type:

numpy.ndarray

ypos

The wavelength-dependent coordinate along the declination direction of the on-sky aperture for all spectra. Expected to be arcseconds from some fiducial coordinate.

Type:

numpy.ndarray

area

The fiducial area subtended by each spectral aperture. The aperture is assumed to be circular.

Type:

numpy.ndarray

prihdr

Primary header for the row-stacked spectra. If not provided on instantiation, set to an empty astropy.io.fits.Header.

Type:

astropy.io.fits.Header

fluxhdr

Header specifically for the flux array. If not provided on instantiation, set to be a copy of prihdr.

Type:

astropy.io.fits.Header

_cube_dimensions(pixelscale=None, recenter=None, width_buffer=None, redo=False)[source]

Determine the on-sky dimensions of the reconstructed image for all wavelength channels.

When first calling this method, the pixel scale, recentering approach, and buffer must be provided, and they are saved to pixelscale, recenter, width_buffer. The spatial dimensions of the cube are calculated using the data in xpos and ypos, and the results are saved to xs, ys, nx, and ny. This calculation only needs to be done once per instance and settings for pixelscale, recenter, and width_buffer. The calculation is forced to be done again if any of these are provided to the method and different from the saved attributes or if redo is True.

Parameters:
  • pixelscale (float, optional) – Desired pixel scale in arcsec.

  • recenter (bool, optional) – Flag to recenter the coordinate system.

  • width_buffer (int, optional) – Number of pixels to use as buffer for the image reconstruction.

  • redo (bool, optional) – Force the recalculation of the cube dimensions if they are already defined.

_formal_covariance_matrix(csr=False, quiet=False)[source]

Return the formal covariance matrix.

The formal covariance matrix is:

\[{\mathbf C} = {\mathbf T} \times {\mathbf \Sigma} \times {\mathbf T}^{\rm T},\]

where \({\mathbf \Sigma}\) is the covariance matrix for the row-stacked spectra for the specified wavelength channel.

This method requires that rect_T not be None, it assumes that the current transfer matrix is correct for this channel, and it is largely a wrapper for mangadap.util.covariance.Covariance.from_matrix_multiplication().

Parameters:
Returns:

mangadap.util.covariance.Covariance, scipy.sparse.csr_matrix: The covariance matrix for the designated wavelength channel. The return type depends on csr.

_init_rectification_parameters(channel=None, pixelscale=None, rlim=None, sigma=None, recenter=None, width_buffer=None)[source]

Perform common preparation of rectification parameters.

Parameters:
  • channel (int, optional) – Index of the spectral channel for which to calculate the transfer matrix.

  • pixelscale (float, optional) – Desired pixel (spaxel) scale in arcsec.

  • rlim (float, optional) – The limiting radius of the image reconstruction (Gaussian) kernel in arcseconds.

  • sigma (float, optional) – The sigma of the image reconstruction (Gaussian) kernel in arcseconds.

  • recenter (bool, optional) – Flag to recenter the coordinate system.

  • width_buffer (int, optional) – Number of pixels to use as buffer for the image reconstruction

_reinit_rectification_parameters(channel=None, pixelscale=None, rlim=None, sigma=None, recenter=None, width_buffer=None)[source]

Determine if the rectification paraemters need to be reinitialized.

binned_on_sky_area(bin_indx, x=None, y=None, **kwargs)[source]

Compute the on-sky area of a set of binned spectra.

The method attempts to calculate the overlapping area assuming the spectra were obtained with a set of circular fibers with the provided area (see area). If not provided in the method call (see x, y), the fiber fiducial coordinates used for the calculation are calculated using mean_sky_coordinates().

The overlapping area is calculated using the shapely python package. If that package is not available, the returned area of each bin is simply the sum of the area of each fiber (does not account for any overlaps).

Parameters:
  • bin_indx (array-like) – An array with size \(N_{\rm spec}\) that gives which spaxels were included in each bin. Valid bins have indices of \(\geq 0\).

  • x (array-like, optional) – On-sky \(x\) coordinate. Default is to calculate \(x\) and \(y\) using mean_sky_coordinates().

  • y (array-like, optional) – On-sky \(y\) coordinate. Default is to calculate \(x\) and \(y\) using mean_sky_coordinates().

  • **kwargs – Arguments passed directly to mean_sky_coordinates() for the determination of the sky coordinates, if they aren’t provided directly.

Returns:

Two numpy.ndarray objects are returned. The first has the unique (non-negative) bin indices, and the second provides the on-sky area of that bin.

Return type:

tuple

central_wavelength(waverange=None, response_func=None, per_pixel=True, flag=None, fluxwgt=False)[source]

Determine the mean central wavelength for all spectra under various conditions.

The wavelength channel is set to be the center of the bandpass selected by waverange, weighted by the response function and the flux (if requested/provided). The mask is also incorporated in these calculations. By default (i.e., no wavelength limits or weighting) the wavelength channel is just the central wavelength of the full spectral range. (Note that if the spectra are binned logarithmically, this isn’t necessarily the central wavelength channel.)

Parameters:
  • waverange (array-like, optional) – Starting and ending wavelength over which to calculate the statistics. Default is to use the full wavelength range.

  • response_func (array-like, optional) – A two-column array with the wavelength and transmission of a broad-band response function to use as a weighting function for the calculation.

  • per_pixel (bool, optional) – When providing a response function, base the calculation on per pixel measurements, instead of per angstrom. Set to False for a per-angstrom calculation.

  • flag (str, list, optional) – One or more flag names that are considered when deciding if a pixel should be masked. The names must be a valid bit name as defined by bitmask. If bitmask is None, these are ignored.

  • fluxwgt (bool, optional) – Flag to weight by the flux when determining the mean coordinates.

Returns:

The mean central wavelength of all spectra.

Return type:

float

copy_to_array(attr='flux', waverange=None, nbins=None, select_bins=None, missing_bins=None, unique_bins=None)[source]

Return a copy of the selected data array with a limited spectral range and/or a selected set of spectra.

The array size is always \(N_{\rm spec} \times N_{\rm wavelength}\). See copy_to_masked_array() for argument descriptions.

Warning

Any masking is ignored in this function call. To incorporate the mask use copy_to_masked_array().

Returns:

A 2D array with a copy of the data from the selected attribute.

Return type:

numpy.ndarray

copy_to_masked_array(attr='flux', use_mask=True, flag=None, waverange=None, nbins=None, select_bins=None, missing_bins=None, unique_bins=None)[source]

Return a copy of the selected data array with selected pixels, wavelength ranges, and/or full spectra masked.

This is functionally identical to copy_to_array(), except the output format is a numpy.ma.MaskedArray. The pixels that are considered to be masked can be specified using the flag option. The array size is always \(N_{\rm spec} \times N_{\rm wavelength}\).

Parameters:
  • attr (str, optional) – The attribute for the returned array. Can be ‘flux’, ‘ivar’, ‘sres’, ‘mask’, ‘xpos’, or ‘ypos’. For ‘mask’, you’re likely better off using copy_to_array(). Strings are always set to be lower-case, so capitalization shouldn’t matter.

  • use_mask (bool, optional) – Use the internal mask to mask the data. This is largely here to allow for copy_to_array() to wrap this function while not applying the internal mask.

  • waverange (array-like, optional) – Two-element array with the first and last wavelength to include in the computation. Default is to use the full wavelength range.

  • flag (str, list, optional) – One or more bitmask flags used to select bad pixels. The names must be a valid bit name as defined by bitmask (see mangadap.util.bitmask.BitMask). If not provided, any non-zero mask bit is omitted.

  • nbins (int, optional) – The total number of defined bins. Default is to find the maximum number in the unique_bins list.

  • select_bins (array-like, optional) – A boolean array selecting the spectra to return.

  • missing_bins (list, optional) – A list of bin numbers that are missing from the output selection and that will be replaced with masked spectra. If specified, must also provide unique_bins.

  • unique_bins (array-like, optional) – The indices of the bins that have valid spectra. The length of this array should match the select_bins array.

Returns:

A 2D array with a copy of the data from the selected extension, masked where mask is either True or with the selected bitmask flags.

Return type:

numpy.ma.MaskedArray

Raises:
  • AttributeError – Raised if attr is not a valid attribute.

  • ValueError – Raised if the unique bins are not specified, but the missing bins are, or if the number of unique bins does not match the number of spectra in the datacube.

covariance_cube(channels=None, pixelscale=None, rlim=None, sigma=None, recenter=False, width_buffer=10, csr=False, quiet=False, rej_flag=None)[source]

Return the covariance matrices for many wavelength channels.

This is primarily a wrapper for covariance_matrix() that iterates through one or more wavelength channels and constructs a 3D mangadap.util.covariance.Covariance object.

Parameters:
  • channels (int, array-like, optional) – Indices of the spectral channels for which to calculate the covariance matrix. If None, the covariance matrix is calculated for all channels.

  • pixelscale (float, optional) – Desired pixel (spaxel) scale in arcsec. If None, pixelscale will be used, but fault if that value is also None.

  • rlim (float, optional) – The limiting radius of the image reconstruction (Gaussian) kernel in arcseconds. If None, rect_rlim will be used, but fault if that value is also None.

  • sigma (float, optional) – The sigma of the image reconstruction (Gaussian) kernel in arcseconds. If None, rect_sigma will be used, but fault if that value is also None.

  • recenter (bool, optional) – Flag to recenter the coordinate system.

  • width_buffer (int, optional) – Number of pixels to use as buffer for the image reconstruction

  • csr (bool, optional) – Instead of returning a mangadap.util.covariance.Covariance object, return the covariance matrix as a scipy.sparse.csr_matrix object. Primarily used by covariance_cube() for collating the covariance matrix of each wavelength channel before combining them into a single mangadap.util.covariance.Covariance object.

  • quiet (bool, optional) – Suppress terminal output

  • rej_flag (str, optional) – Bit name to ignore in rectification. Ignored if bitmask is None. If None, bitmasks are ignored. Use 'any' to mask pixels with any flagged bits.

Returns:

mangadap.util.covariance.Covariance, numpy.ndarray: If csr is True, the returned object is an numpy.ndarray of scipy.sparse.csr_matrix types.

covariance_matrix(channel, pixelscale=None, rlim=None, sigma=None, recenter=False, width_buffer=10, csr=False, quiet=False, rej_flag=None)[source]

Return the covariance matrix for the specified wavelength channel.

For a regrided cube image with \(N_x\times N_y\) pixels, the covariance matrix has a size that is \((N_x N_y)\times (N_x N_y)\); however, the majority of these pixels will be zero. Therefore, the covariance matrix is stored as a sparse matrix and interfaced with using the mangadap.util.covariance.Covariance object class.

The value of the covariance matrix at pixel \((i,j)\) is the covariance between pixels \((n_{x,0},n_{y,0})\) and \((n_{x,1},n_{y,1})\) at the specified wavelength channel of the reconstructed datacube image, where

\[\begin{split}n_{x,0} &= \lfloor i / N_y \rfloor \\ n_{y,0} &= i - n_{x,0} N_y \\ n_{x,1} &= \lfloor j / N_y \rfloor \\ n_{y,1} &= j - n_{x,1} N_y\end{split}\]

and \(\lfloor m\rfloor\) is the “floor” of \(m\). The diagonal of the covariance matrix (\(i=j\)) provides the variance.

Warning

Because the covariance is calculated for a higher dimensional array, it’s important to note how the pixels in the image map to the relevant covariance array locations. See the equations above. Also see mangadap.util.covariance.Covariance.transpose_raw_shape().

The covariance matrix provided by this method is always the formal covariance matrix defined as

\[{\mathbf C} = {\mathbf T} \times {\mathbf \Sigma} \times {\mathbf T}^{\rm T},\]

where \({\mathbf \Sigma}\) is the covariance matrix for the row-stacked spectra for the specified wavelength channel, which is taken to be a diagonal matrix with the flux variances along the diagonal.

Parameters:
  • channel (int) – Index of the spectral channel for which to calculate the covariance matrix.

  • pixelscale (float, optional) – Desired pixel (spaxel) scale in arcsec. If None, pixelscale will be used, but fault if that value is also None.

  • rlim (float, optional) – The limiting radius of the image reconstruction (Gaussian) kernel in arcseconds. If None, rect_rlim will be used, but fault if that value is also None.

  • sigma (float, optional) – The sigma of the image reconstruction (Gaussian) kernel in arcseconds. If None, rect_sigma will be used, but fault if that value is also None.

  • recenter (bool, optional) – Flag to recenter the coordinate system.

  • width_buffer (int, optional) – Number of pixels to use as buffer for the image reconstruction

  • csr (bool, optional) – Instead of returning a mangadap.util.covariance.Covariance object, return the covariance matrix as a scipy.sparse.csr_matrix object. Primarily used by covariance_cube() for collating the covariance matrix of each wavelength channel before combining them into a single mangadap.util.covariance.Covariance object.

  • quiet (bool, optional) – Suppress terminal output

  • rej_flag (str, optional) – Bit name to ignore in rectification. Ignored if bitmask is None. If None, bitmasks are ignored. Use 'any' to mask pixels with any flagged bits.

Returns:

mangadap.util.covariance.Covariance, scipy.sparse.csr_matrix: The covariance matrix for the designated wavelength channel. The return type depends on csr.

flux_stats(waverange=None, response_func=None, per_pixel=True, flag=None)[source]

Compute the mean flux, propagated error in the mean flux, and mean S/N over the specified wavelength range.

If the wavelength range is not specified, the quantities are calculated over the full spectral range.

Parameters:
  • waverange (array-like, optional) – Starting and ending wavelength over which to calculate the statistics. Default is to use the full wavelength range.

  • response_func (array-like, optional) – A two-column array with the wavelength and transmission of a broad-band response function to use as a weighting function for the calculation.

  • per_pixel (bool, optional) – When providing a response function, continue to calculate the statistics per pixel. Set to False for a per-angstrom calculation.

  • flag (str, list, optional) – One or more flag names that are considered when deciding if a pixel should be masked. The names must be a valid bit name as defined by bitmask.

Returns:

Three objects are returned: the mean flux, the propagated variance in the mean flux, and the mean S/N.

Return type:

numpy.ndarray

Raises:

ValueError – Raised of a provided wavelength range object does not have two elements or if a response function is provided and has an incorrect shape.

instrumental_dispersion_plane(channel, dispersion_factor=None, pixelscale=None, rlim=None, sigma=None, recenter=None, width_buffer=None, quiet=False, rej_flag=None)[source]

Return the instrumental dispersion for the rectified wavelength plane.

Method requires that sres is not None.

The method first calculates the rectification transfer matrix, \({\mathbf T}\), using regrid_transfer_matrix(). The transfer matrix is then used to construct the rectified wavelength plane image, \({\mathbf I}\), by computing \({\mathbf T} \times {\mathbf F} = {\mathbf I}\), where \({\mathbf F}\) is the vector of the fiber fluxes. Under the assumption that the line-spread-function (LSF) is Gaussian, we determine the instrumental dispersion for the data in the wavelength channel of the reconstructed CUBE by calculating the second moment of the weighted sum of Gaussians of the appropriate dispersion. Assuming all the Gaussians have the normalized form:

\[\mathcal{N}(x|\mu=0,\sigma) = \frac{1}{\sqrt{2\pi}\sigma} \exp\left(-\frac{x^2}{2\sigma^2}\right),\]

the combined instrumental dispersion becomes

\[\sigma_{\rm inst,j}^2 = \frac{\sum_i T_{ij} \sigma^2_i}{\sum_i T_{ij}},\]

where \(T_{ij}\) are the elements of the transfer matrix. In terms of matrix multiplications, this can be written as

\[{\mathbf S} = \frac{ {\mathbf T} \times {\mathbf V} }{ {\mathbf T_c} },\]

where \({\mathbf T_c} = {\mathbf T_c} \times {\mathbf 1}\) is the vector with the sum \(\sum_j T_{ij}\), \({\mathbf V}\) is the instrumental variance for all fibers at the designated wavelength plane, and \({\mathbf S}\) is the variance for all the spaxels in the reconstructed wavelength image; the division by \({\mathbf T_c}\) is element-wise.

The returned matrix is the element-wise square-root of \({\mathbf S}\), rearranged into a 2D array of size nx by ny.

Parameters:
  • channel (int) – Index of the spectral channel for which to calculate the covariance matrix.

  • dispersion_factor (float, optional) – Artificially change the resolution by multiplying the instrumental LSF by this factor before calculating the reconstructed dispersion. If None, factor is unity.

  • pixelscale (float, optional) – Desired pixel (spaxel) scale in arcsec. If None, pixelscale will be used, but fault if that value is also None.

  • rlim (float, optional) – The limiting radius of the image reconstruction (Gaussian) kernel in arcseconds. If None, rect_rlim will be used, but fault if that value is also None.

  • sigma (float, optional) – The sigma of the image reconstruction (Gaussian) kernel in arcseconds. If None, rect_sigma will be used, but fault if that value is also None.

  • recenter (bool, optional) – Flag to recenter the coordinate system.

  • width_buffer (int, optional) – Number of pixels to use as buffer for the image reconstruction

  • quiet (bool, optional) – Suppress terminal output

  • rej_flag (str, optional) – Bit name to ignore in rectification. Ignored if bitmask is None. If None, bitmasks are ignored. Use 'any' to mask pixels with any flagged bits.

Returns:

The reconstructed LSF dispersion for the specified wavelength channel.

Return type:

numpy.ndarray

interpolate_to_match(func, fill_value=0.0)[source]

Interpolate a function to match the datacube wavelength sampling.

Parameters:
  • func (numpy.ndarray) – Function to linear interpolate. Shape is \((N_{\rm wave}, 2)\), where the first column has the wavelength and the second has the function to interpolate. If None, simply returns a unity vector of the correct length.

  • fill_value (float, optional) – The value to use for spectral regions not sampled by the provided function.

Returns:

The function interpolated (not resampled) to match wave. Any spectral regions not within the wavelength range of the provided function is set to 0.

Return type:

numpy.ndarray

mean_sky_coordinates(waverange=None, response_func=None, per_pixel=True, offset=None, flag=None, fluxwgt=False)[source]

Compute the weighted or unweighted mean sky coordinates for each spectrum.

Method cannot be performed if coordinate arrays are not available.

Warning

Flux-weighting the coordinates can produce spurious results in low-flux regimes.

Parameters:
  • waverange (array-like, optional) – Two-element array with the first and last wavelength to include in the computation. Default is to use the full wavelength range.

  • response_func (array-like, optional) – A two-column array with the wavelength and transmission of a broad-band response function to use for the calculation. If None, no wavelength-dependent weighting function is used.

  • per_pixel (bool, optional) – When providing a response function, continue to calculate the statistics per pixel, as opposed to per angstrom.

  • offset (tuple, optional) – A two-tuple with an x,y offset to apply.

  • flag (str, list, optional) – One or more flag names that are considered when deciding if a pixel should be masked. The names must be a valid bit name as defined by bitmask.

  • fluxwgt (bool, optional) – Flag to weight by the flux when determining the mean coordinates.

Returns:

Two numpy.ndarray objects with the fiducial x and y on-sky positions of each spectrum in arcseconds relative to a given center. The shape of the two arrays is \((N_{\rm spec},)\).

Return type:

tuple

property nspec

Number of spectra in the datacube.

rectification_transfer_matrix(channel, pixelscale, rlim, sigma, recenter=False, width_buffer=10, quiet=False, rej_flag=None)[source]

Calculate the transfer matrix used to produce a reconstructed image of the fiber data at the specified wavelength channel using Shepard’s method.

This is done directly using the available on-sky x and y coordinates of each fiber as a function of wavelength taken from xpos and ypos attributes.

Todo

  • Give more detail on the pixels at which the radius is calculated.

  • Describe the algorithm in more depth.

Parameters:
  • channel (int) – Index of the spectral channel for which to calculate the transfer matrix.

  • pixelscale (float) – Desired pixel (spaxel) scale in arcsec.

  • rlim (float) – The limiting radius of the image reconstruction (Gaussian) kernel in arcseconds.

  • sigma (float) – The sigma of the image reconstruction (Gaussian) kernel in arcseconds.

  • recenter (bool, optional) – Flag to recenter the coordinate system.

  • width_buffer (int, optional) – Number of pixels to use as buffer for the image reconstruction

  • quiet (bool, optional) – Suppress terminal output

  • rej_flag (str, optional) – Bit name to ignore in rectification. Ignored if bitmask is None. If None, bitmasks are ignored. Use 'any' to mask pixels with any flagged bits.

Returns:

Transfer matrix \({\mathbf T}\)

Return type:

scipy.sparse.csr_matrix

rectify_wavelength_plane(channel, pixelscale=None, rlim=None, sigma=None, recenter=False, width_buffer=10, quiet=False, rej_flag=None, return_ivar=False, return_covar=False)[source]

Return a rectified image for the specified wavelength channel using Shepard’s method.

First, the rectification transfer matrix, \({\mathbf T}\), is calculated using rectification_transfer_matrix(). Then the wavelength image, \({\mathbf I}\) is calculated:

\[{\mathbf T} \times {\mathbf F} = {\mathbf I}\]

where \({\mathbf F}\) is the vector of fluxes in the selected wavelength channel for all the aperture measurements. The calculation of \({\mathbf I}\) produces a vector, but this is reshaped into a 2D array of shape \((N_x, N_y)\) on output.

Parameters:
  • channel (int) – Index of the spectral channel for which to calculate the transfer matrix.

  • pixelscale (float, optional) – Desired pixel (spaxel) scale in arcsec. If None, pixelscale will be used, but fault if that value is also None.

  • rlim (float, optional) – The limiting radius of the image reconstruction (Gaussian) kernel in arcseconds. If None, rect_rlim will be used, but fault if that value is also None.

  • sigma (float, optional) – The sigma of the image reconstruction (Gaussian) kernel in arcseconds. If None, rect_sigma will be used, but fault if that value is also None.

  • recenter (bool, optional) – Flag to recenter the coordinate system.

  • width_buffer (int, optional) – Number of pixels to use as buffer for the image reconstruction

  • quiet (bool, optional) – Suppress terminal output

  • rej_flag (str, optional) – Bit name to ignore in rectification. Ignored if bitmask is None. If None, bitmasks are ignored. Use 'any' to mask pixels with any flagged bits.

  • return_ivar (bool, optional) – Return the nominal inverse variance image (i.e. the diagonal of the covariance matrix).

  • return_covar (bool, optional) – Return the covariance matrix. Note that the covariance matrix is calculated if either return_ivar or return_covar are True. If return_covar is True, return_ivar is ignored (because it’s automatically true).

Returns:

If both return_ivar and return_covar are False, only the reconstructed image is returned. If return_covar is True, a mangadap.util.covariance.Covariance object is also returned. If return_covar is False and return_ivar is True, the second returned object is a numpy.ndarray with the propagated inverse variance image.

Return type:

numpy.ndarray, tuple