pypeit.core.datacube module
Module containing routines used by 3D datacubes.
- class pypeit.core.datacube.DataCube(flux, sig, bpm, PYP_SPEC, blaze_wave, blaze_spec, sensfunc=None, fluxed=None)[source]
Bases:
DataContainer
DataContainer to hold the products of a datacube
The datamodel attributes are:
Version: 1.1.0
Attribute
Type
Array Type
Description
PYP_SPEC
str
PypeIt: Spectrograph name
blaze_spec
The spectral blaze function
blaze_wave
Wavelength array of the spectral blaze function
bpm
Bad pixel mask of the datacube (0=good, 1=bad)
flux
Flux datacube in units of counts/s/Ang/arcsec^2 or 10^-17 erg/s/cm^2/Ang/arcsec^2
fluxed
bool
Boolean indicating if the datacube is fluxed.
sensfunc
Sensitivity function 10^-17 erg/(counts/cm^2)
sig
Error datacube (matches units of flux)
- Parameters:
flux (numpy.ndarray) – The science datacube (nwave, nspaxel_y, nspaxel_x)
sig (numpy.ndarray) – The error datacube (nwave, nspaxel_y, nspaxel_x)
bpm (numpy.ndarray) – The bad pixel mask of the datacube (nwave, nspaxel_y, nspaxel_x). True values indicate a bad pixel
blaze_wave (numpy.ndarray) – Wavelength array of the spectral blaze function
blaze_spec (numpy.ndarray) – The spectral blaze function
sensfunc (numpy.ndarray, None) – Sensitivity function (nwave,). Only saved if the data are fluxed.
PYP_SPEC (str) – Name of the PypeIt Spectrograph
fluxed (bool) – If the cube has been flux calibrated, this will be set to “True”
- head0
Primary header
- Type:
- spectrograph
Build from PYP_SPEC
- Type:
- _bundle()[source]
Over-write default _bundle() method to separate the DetectorContainer into its own HDU
- Returns:
A list of dictionaries, each list element is written to its own fits extension. See the description above.
- Return type:
- datamodel = {'PYP_SPEC': {'descr': 'PypeIt: Spectrograph name', 'otype': <class 'str'>}, 'blaze_spec': {'atype': <class 'numpy.floating'>, 'descr': 'The spectral blaze function', 'otype': <class 'numpy.ndarray'>}, 'blaze_wave': {'atype': <class 'numpy.floating'>, 'descr': 'Wavelength array of the spectral blaze function', 'otype': <class 'numpy.ndarray'>}, 'bpm': {'atype': <class 'numpy.uint8'>, 'descr': 'Bad pixel mask of the datacube (0=good, 1=bad)', 'otype': <class 'numpy.ndarray'>}, 'flux': {'atype': <class 'numpy.floating'>, 'descr': 'Flux datacube in units of counts/s/Ang/arcsec^2 or 10^-17 erg/s/cm^2/Ang/arcsec^2', 'otype': <class 'numpy.ndarray'>}, 'fluxed': {'descr': 'Boolean indicating if the datacube is fluxed.', 'otype': <class 'bool'>}, 'sensfunc': {'atype': <class 'numpy.floating'>, 'descr': 'Sensitivity function 10^-17 erg/(counts/cm^2)', 'otype': <class 'numpy.ndarray'>}, 'sig': {'atype': <class 'numpy.floating'>, 'descr': 'Error datacube (matches units of flux)', 'otype': <class 'numpy.ndarray'>}}
Provides the class data model. This is undefined in the abstract class and should be overwritten in the derived classes.
The format of the
datamodel
needed for each implementation of aDataContainer
derived class is as follows.The datamodel itself is a class attribute (i.e., it is a member of the class, not just of an instance of the class). The datamodel is a dictionary of dictionaries: Each key of the datamodel dictionary provides the name of a given datamodel element, and the associated item (dictionary) for the datamodel element provides the type and description information for that datamodel element. For each datamodel element, the dictionary item must provide:
otype
: This is the type of the object for this datamodel item. E.g., for a float or a numpy.ndarray, you would setotype=float
andotype=np.ndarray
, respectively.descr
: This provides a text description of the datamodel element. This is used to construct the datamodel tables in the pypeit documentation.
If the object type is a numpy.ndarray, you should also provide the
atype
keyword that sets the type of the data contained within the array. E.g., for a floating point array containing an image, your datamodel could be simply:datamodel = {'image' : dict(otype=np.ndarray, atype=float, descr='My image')}
More advanced examples are given in the top-level module documentation.
Currently,
datamodel
components are restricted to haveotype
that aretuple
,int
,float
,numpy.integer
,numpy.floating
, numpy.ndarray, or astropy.table.Table objects. E.g.,datamodel
values forotype
cannot bedict
.
- classmethod from_file(ifile)[source]
Over-load
from_file()
to deal with the header- Parameters:
ifile (str) – Filename holding the object
- internals = ['head0', 'filename', 'spectrograph', 'spect_meta']
A list of strings identifying a set of internal attributes that are not part of the datamodel.
- property ivar
- to_file(ofile, primary_hdr=None, hdr=None, **kwargs)[source]
Over-load
to_file()
to deal with the header- Parameters:
ofile (
str
) – Filenameprimary_hdr (astropy.io.fits.Header, optional) – Base primary header. Updated with new subheader keywords. If None, initialized using
initialize_header()
.wcs (astropy.io.fits.Header, optional) – The World Coordinate System, represented by a fits header
kwargs (dict) – Keywords passed directly to parent
to_file
function.
- version = '1.1.0'
Provides the string representation of the class version.
This is currently put to minimal use so far, but will used for I/O verification in the future.
Each derived class should provide a version to guard against data model changes during development.
- property wcs
- pypeit.core.datacube.coadd_cube(files, opts, spectrograph=None, parset=None, overwrite=False)[source]
Main routine to coadd spec2D files into a 3D datacube
- Parameters:
files (
list
) – List of all spec2D filesopts (
dict
) – coadd2d options associated with each spec2d filespectrograph (
str
,Spectrograph
, optional) – The name or instance of the spectrograph used to obtain the data. If None, this is pulled from the file header.parset (
PypeItPar
, optional) – An instance of the parameter set. If None, assumes that detector 1 is the one reduced and uses the default reduction parameters for the spectrograph (seedefault_pypeit_par()
for the relevant spectrograph class).overwrite (
bool
, optional) – Overwrite the output file, if it exists?
- pypeit.core.datacube.compute_weights(all_ra, all_dec, all_wave, all_sci, all_ivar, all_idx, whitelight_img, dspat, dwv, sn_smooth_npix=None, relative_weights=False)[source]
Calculate wavelength dependent optimal weights. The weighting is currently based on a relative \((S/N)^2\) at each wavelength
- Parameters:
all_ra (numpy.ndarray) – 1D flattened array containing the RA values of each pixel from all spec2d files
all_dec (numpy.ndarray) – 1D flattened array containing the DEC values of each pixel from all spec2d files
all_wave (numpy.ndarray) – 1D flattened array containing the wavelength values of each pixel from all spec2d files
all_sci (numpy.ndarray) – 1D flattened array containing the counts of each pixel from all spec2d files
all_ivar (numpy.ndarray) – 1D flattened array containing the inverse variance of each pixel from all spec2d files
all_idx (numpy.ndarray) – 1D flattened array containing an integer identifier indicating which spec2d file each pixel originates from. For example, a 0 would indicate that a pixel originates from the first spec2d frame listed in the input file. a 1 would indicate that this pixel originates from the second spec2d file, and so forth.
whitelight_img (numpy.ndarray) – A 2D array containing a whitelight image, that was created with the input
all_
arrays.dspat (float) – The size of each spaxel on the sky (in degrees)
dwv (float) – The size of each wavelength pixel (in Angstroms)
sn_smooth_npix (float, optional) – Number of pixels used for determining smoothly varying S/N ratio weights. This is currently not required, since a relative weighting scheme with a polynomial fit is used to calculate the S/N weights.
relative_weights (bool, optional) – Calculate weights by fitting to the ratio of spectra?
- Returns:
a 1D array the same size as all_sci, containing relative wavelength dependent weights of each input pixel.
- Return type:
- pypeit.core.datacube.correct_dar(wave_arr, coord, obstime, location, pressure, temperature, rel_humidity, wave_ref=None, numgrid=10)[source]
Apply a differental atmospheric refraction correction to the input ra/dec.
This implementation is based on ERFA, which is called through astropy.
Todo
There’s probably going to be issues when the RA angle is either side of RA=0.
- Parameters:
wave_arr (numpy.ndarray) – wavelengths to obtain ra and dec offsets
coord (astropy.coordinates.SkyCoord) – ra, dec positions at the centre of the field
obstime (astropy.time.Time) – time at the midpoint of observation
location (astropy.coordinates.EarthLocation) – observatory location
pressure (
float
) – Outside pressure at locationtemperature (
float
) – Outside ambient air temperature at locationrel_humidity (
float
) – Outside relative humidity at location. This should be between 0 to 1.wave_ref (
float
) – Reference wavelength (The DAR correction will be performed relative to this wavelength)numgrid (
int
) – Number of grid points to evaluate the DAR correction.
- Returns:
ra_diff (numpy.ndarray) – Relative RA shift at each wavelength given by
wave_arr
dec_diff (numpy.ndarray) – Relative DEC shift at each wavelength given by
wave_arr
- pypeit.core.datacube.correct_grating_shift(wave_eval, wave_curr, spl_curr, wave_ref, spl_ref, order=2)[source]
Using spline representations of the blaze profile, calculate the grating correction that should be applied to the current spectrum (suffix
curr
) relative to the reference spectrum (suffixref
). The grating correction is then evaluated at the wavelength array given bywave_eval
.- Parameters:
wave_eval (numpy.ndarray) – Wavelength array to evaluate the grating correction
wave_curr (numpy.ndarray) – Wavelength array used to construct spl_curr
spl_curr (scipy.interpolate.interp1d) – Spline representation of the current blaze function (based on the illumflat).
wave_ref (numpy.ndarray) – Wavelength array used to construct spl_ref
spl_ref (scipy.interpolate.interp1d) – Spline representation of the reference blaze function (based on the illumflat).
order (int) – Polynomial order used to fit the grating correction.
- Returns:
The grating correction to apply
- Return type:
- pypeit.core.datacube.create_wcs(cubepar, all_ra, all_dec, all_wave, dspat, dwv, collapse=False, equinox=2000.0, specname='PYP_SPEC')[source]
Create a WCS and the expected edges of the voxels, based on user-specified parameters or the extremities of the data.
- Parameters:
cubepar (
CubePar
) – An instance of the CubePar parameter set, contained parameters of the datacube reductionall_ra (numpy.ndarray) – 1D flattened array containing the RA values of each pixel from all spec2d files
all_dec (numpy.ndarray) – 1D flattened array containing the DEC values of each pixel from all spec2d files
all_wave (numpy.ndarray) – 1D flattened array containing the wavelength values of each pixel from all spec2d files
dspat (float) – Spatial size of each square voxel (in arcsec). The default is to use the values in cubepar.
dwv (float) – Linear wavelength step of each voxel (in Angstroms)
collapse (bool, optional) – If True, the spectral dimension will be collapsed to a single channel (primarily for white light images)
equinox (float, optional) – Equinox of the WCS
specname (str, optional) – Name of the spectrograph
- Returns:
cubewcs (astropy.wcs.WCS) – astropy WCS to be used for the combined cube
voxedges (tuple) – A three element tuple containing the bin edges in the x, y (spatial) and z (wavelength) dimensions
reference_image (numpy.ndarray) – The reference image to be used for the cross-correlation. Can be None.
- pypeit.core.datacube.dar_fitfunc(radec, coord_ra, coord_dec, datfit, wave, obstime, location, pressure, temperature, rel_humidity)[source]
Generates a fitting function to calculate the offset due to differential atmospheric refraction
- Parameters:
radec (tuple) – A tuple containing two floats representing the shift in ra and dec due to DAR.
coord_ra (float) – RA in degrees
coord_dec (float) – Dec in degrees
datfit (numpy.ndarray) – The RA and DEC that the model needs to match
wave (float) – Wavelength to calculate the DAR
location (astropy.coordinates.EarthLocation) – observatory location
pressure (float) – Outside pressure at location
temperature (float) – Outside ambient air temperature at location
rel_humidity (float) – Outside relative humidity at location. This should be between 0 to 1.
- Returns:
chi-squared difference between datfit and model
- Return type:
- pypeit.core.datacube.extract_standard_spec(stdcube, subpixel=20, method='boxcar')[source]
Extract a spectrum of a standard star from a datacube
- Parameters:
std_cube (astropy.io.fits.HDUList) – An HDU list of fits files
subpixel (int) – Number of pixels to subpixelate spectrum when creating mask
method (str) – Method used to extract standard star spectrum. Currently, only ‘boxcar’ is supported
- Returns:
wave (numpy.ndarray) – Wavelength of the star.
Nlam_star (numpy.ndarray) – counts/second/Angstrom
Nlam_ivar_star (numpy.ndarray) – inverse variance of Nlam_star
gpm_star (numpy.ndarray) – good pixel mask for Nlam_star
- pypeit.core.datacube.fitGaussian2D(image, norm=False)[source]
Fit a 2D Gaussian to an input image. It is recommended that the input image is scaled to a maximum value that is ~1, so that all fit parameters are of the same order of magnitude. Set norm=True if you do not care about the amplitude or integrated flux. Otherwise, make sure you scale the image by a known value prior to passing it into this function.
- Parameters:
image (numpy.ndarray) – A 2D input image
norm (bool, optional) – If True, the input image will be normalised to the maximum value of the input image.
- Returns:
popt (numpy.ndarray) – The optimum parameters of the Gaussian in the following order: Integrated flux, x center, y center, sigma_x, sigma_y, theta, offset. See
gaussian2D()
for a more detailed description of the model.pcov (numpy.ndarray) – Corresponding covariance matrix
- pypeit.core.datacube.gaussian2D(tup, intflux, xo, yo, sigma_x, sigma_y, theta, offset)[source]
Fit a 2D Gaussian function to an image.
- Parameters:
tup (
tuple
) – A two element tuple containing the x and y coordinates of each pixel in the imageintflux (float) – The Integrated flux of the 2D Gaussian
xo (float) – The centre of the Gaussian along the x-coordinate when z=0
yo (float) – The centre of the Gaussian along the y-coordinate when z=0
sigma_x (float) – The standard deviation in the x-direction
sigma_y (float) – The standard deviation in the y-direction
theta (float) – The orientation angle of the 2D Gaussian
offset (float) – Constant offset
- Returns:
The 2D Gaussian evaluated at the coordinate (x, y)
- Return type:
- pypeit.core.datacube.gaussian2D_cube(tup, intflux, xo, yo, dxdz, dydz, sigma_x, sigma_y, theta, offset)[source]
Fit a 2D Gaussian function to a datacube. This function assumes that each wavelength slice of the datacube is well-fit by a 2D Gaussian. The centre of the Gaussian is allowed to vary linearly as a function of wavelength.
Note
The integrated flux does not vary with wavelength.
- Parameters:
tup (
tuple
) – A three element tuple containing the x, y, and z locations of each pixel in the cubeintflux (float) – The Integrated flux of the Gaussian
xo (float) – The centre of the Gaussian along the x-coordinate when z=0
yo (float) – The centre of the Gaussian along the y-coordinate when z=0
dxdz (float) – The change of xo with increasing z
dydz (float) – The change of yo with increasing z
sigma_x (float) – The standard deviation in the x-direction
sigma_y (float) – The standard deviation in the y-direction
theta (float) – The orientation angle of the 2D Gaussian
offset (float) – Constant offset
- Returns:
The 2D Gaussian evaluated at the coordinate (x, y, z)
- Return type:
- pypeit.core.datacube.generate_WCS(crval, cdelt, equinox=2000.0, name='PYP_SPEC')[source]
Generate a WCS that will cover all input spec2D files
- Parameters:
- Returns:
astropy WCS to be used for the combined cube
- Return type:
- pypeit.core.datacube.generate_cube_subpixel(outfile, output_wcs, all_ra, all_dec, all_wave, all_sci, all_ivar, all_wghts, all_spatpos, all_specpos, all_spatid, tilts, slits, astrom_trans, bins, all_idx=None, spec_subpixel=10, spat_subpixel=10, overwrite=False, blaze_wave=None, blaze_spec=None, fluxcal=False, sensfunc=None, whitelight_range=None, specname='PYP_SPEC', debug=False)[source]
Save a datacube using the subpixel algorithm. Refer to the subpixellate() docstring for further details about this algorithm
- Parameters:
outfile (str) – Filename to be used to save the datacube
output_wcs (astropy.wcs.WCS) – Output world coordinate system.
all_ra (numpy.ndarray) – 1D flattened array containing the right ascension of each pixel (units = degrees)
all_dec (numpy.ndarray) – 1D flattened array containing the declination of each pixel (units = degrees)
all_wave (numpy.ndarray) – 1D flattened array containing the wavelength of each pixel (units = Angstroms)
all_sci (numpy.ndarray) – 1D flattened array containing the counts of each pixel from all spec2d files
all_ivar (numpy.ndarray) – 1D flattened array containing the inverse variance of each pixel from all spec2d files
all_wghts (numpy.ndarray) – 1D flattened array containing the weights of each pixel to be used in the combination
all_spatpos (numpy.ndarray) – 1D flattened array containing the detector pixel location in the spatial direction
all_specpos (numpy.ndarray) – 1D flattened array containing the detector pixel location in the spectral direction
all_spatid (numpy.ndarray) – 1D flattened array containing the spatid of each pixel
tilts (numpy.ndarray, list) – 2D wavelength tilts frame, or a list of tilt frames (see all_idx)
slits (
SlitTraceSet
, list) – Information stored about the slits, or a list of SlitTraceSet (see all_idx)astrom_trans (
AlignmentSplines
, list) – A Class containing the transformation between detector pixel coordinates and WCS pixel coordinates, or a list of Alignment Splines (see all_idx)bins (tuple) – A 3-tuple (x,y,z) containing the histogram bin edges in x,y spatial and z wavelength coordinates
all_idx (numpy.ndarray, optional) – If tilts, slits, and astrom_trans are lists, this should contain a 1D flattened array, of the same length as all_sci, containing the index the tilts, slits, and astrom_trans lists that corresponds to each pixel. Note that, in this case all of these lists need to be the same length.
spec_subpixel (int, optional) – What is the subpixellation factor in the spectral direction. Higher values give more reliable results, but note that the time required goes as (
spec_subpixel * spat_subpixel
). The default value is 5, which divides each detector pixel into 5 subpixels in the spectral direction.spat_subpixel (int, optional) – What is the subpixellation factor in the spatial direction. Higher values give more reliable results, but note that the time required goes as (
spec_subpixel * spat_subpixel
). The default value is 5, which divides each detector pixel into 5 subpixels in the spatial direction.overwrite (bool, optional) – If True, the output cube will be overwritten.
blaze_wave (numpy.ndarray, optional) – Wavelength array of the spectral blaze function
blaze_spec (numpy.ndarray, optional) – Spectral blaze function
fluxcal (bool, optional) – Are the data flux calibrated? If True, the units are: \({\rm erg/s/cm}^2{\rm /Angstrom/arcsec}^2\) multiplied by the PYPEIT_FLUX_SCALE. Otherwise, the units are: \({\rm counts/s/Angstrom/arcsec}^2\).
sensfunc (numpy.ndarray, None, optional) – Sensitivity function that has been applied to the datacube
whitelight_range (None, list, optional) – A two element list that specifies the minimum and maximum wavelengths (in Angstroms) to use when constructing the white light image (format is: [min_wave, max_wave]). If None, the cube will be collapsed over the full wavelength range. If a list is provided an either element of the list is None, then the minimum/maximum wavelength range of that element will be set by the minimum/maximum wavelength of all_wave.
specname (str, optional) – Name of the spectrograph
debug (bool, optional) – If True, a residuals cube will be output. If the datacube generation is correct, the distribution of pixels in the residual cube with no flux should have mean=0 and std=1.
- pypeit.core.datacube.generate_image_subpixel(image_wcs, all_ra, all_dec, all_wave, all_sci, all_ivar, all_wghts, all_spatpos, all_specpos, all_spatid, tilts, slits, astrom_trans, bins, all_idx=None, spec_subpixel=10, spat_subpixel=10, combine=False)[source]
Generate a white light image from the input pixels
- Parameters:
image_wcs (astropy.wcs.WCS) – World coordinate system to use for the white light images.
all_ra (numpy.ndarray) – 1D flattened array containing the right ascension of each pixel (units = degrees)
all_dec (numpy.ndarray) – 1D flattened array containing the declination of each pixel (units = degrees)
all_wave (numpy.ndarray) – 1D flattened array containing the wavelength of each pixel (units = Angstroms)
all_sci (numpy.ndarray) – 1D flattened array containing the counts of each pixel from all spec2d files
all_ivar (numpy.ndarray) – 1D flattened array containing the inverse variance of each pixel from all spec2d files
all_wghts (numpy.ndarray) – 1D flattened array containing the weights of each pixel to be used in the combination
all_spatpos (numpy.ndarray) – 1D flattened array containing the detector pixel location in the spatial direction
all_specpos (numpy.ndarray) – 1D flattened array containing the detector pixel location in the spectral direction
all_spatid (numpy.ndarray) – 1D flattened array containing the spatid of each pixel
tilts (numpy.ndarray, list) – 2D wavelength tilts frame, or a list of tilt frames (see all_idx)
slits (
SlitTraceSet
, list) – Information stored about the slits, or a list of SlitTraceSet (see all_idx)astrom_trans (
AlignmentSplines
, list) – A Class containing the transformation between detector pixel coordinates and WCS pixel coordinates, or a list of Alignment Splines (see all_idx)bins (tuple) – A 3-tuple (x,y,z) containing the histogram bin edges in x,y spatial and z wavelength coordinates
all_idx (numpy.ndarray, optional) – If tilts, slits, and astrom_trans are lists, this should contain a 1D flattened array, of the same length as all_sci, containing the index the tilts, slits, and astrom_trans lists that corresponds to each pixel. Note that, in this case all of these lists need to be the same length.
spec_subpixel (
int
, optional) – What is the subpixellation factor in the spectral direction. Higher values give more reliable results, but note that the time required goes as (spec_subpixel * spat_subpixel
). The default value is 5, which divides each detector pixel into 5 subpixels in the spectral direction.spat_subpixel (
int
, optional) – What is the subpixellation factor in the spatial direction. Higher values give more reliable results, but note that the time required goes as (spec_subpixel * spat_subpixel
). The default value is 5, which divides each detector pixel into 5 subpixels in the spatial direction.combine (
bool
, optional) – If True, all of the input frames will be combined into a single output. Otherwise, individual images will be generated.
- Returns:
The white light images for all frames
- Return type:
- pypeit.core.datacube.get_output_filename(fil, par_outfile, combine, idx=1)[source]
Get the output filename of a datacube, given the input
- Parameters:
- Returns:
The output filename to use.
- Return type:
- pypeit.core.datacube.get_output_whitelight_filename(outfile)[source]
Given the output filename of a datacube, create an appropriate whitelight fits file name
- pypeit.core.datacube.get_whitelight_pixels(all_wave, min_wl, max_wl)[source]
Determine which pixels are included within the specified wavelength range
- Parameters:
all_wave (numpy.ndarray) – The wavelength of each individual pixel
min_wl (float) – Minimum wavelength to consider
max_wl (float) – Maximum wavelength to consider
- Returns:
A numpy.ndarray object with the indices of all_wave that contain pixels within the requested wavelength range, and a float with the wavelength range (i.e. maximum wavelength - minimum wavelength)
- Return type:
- pypeit.core.datacube.get_whitelight_range(wavemin, wavemax, wl_range)[source]
Get the wavelength range to use for the white light images
- Parameters:
wavemin (float) – Automatically determined minimum wavelength to use for making the white light image.
wavemax (float) – Automatically determined maximum wavelength to use for making the white light image.
wl_range (list) – Two element list containing the user-specified values to manually override the automated values determined by PypeIt.
- Returns:
wlrng – A two element list containing the minimum and maximum wavelength to use for the white light images
- Return type:
- pypeit.core.datacube.load_imageWCS(filename, ext=0)[source]
Load an image and return the image and the associated WCS.
- Parameters:
- Returns:
An numpy.ndarray with the 2D image data and a astropy.wcs.WCS with the image WCS.
- Return type:
- pypeit.core.datacube.make_good_skymask(slitimg, tilts)[source]
Mask the spectral edges of each slit (i.e. the pixels near the ends of the detector in the spectral direction). Some extreme values of the tilts are only sampled with a small fraction of the pixels of the slit width. This leads to a bad extrapolation/determination of the sky model.
- Parameters:
slitimg (numpy.ndarray) – An image of the slit indicating which slit each pixel belongs to
tilts (numpy.ndarray) – Spectral tilts.
- Returns:
A mask of the good sky pixels (True = good)
- Return type:
- pypeit.core.datacube.make_whitelight_fromcube(cube, wave=None, wavemin=None, wavemax=None)[source]
Generate a white light image using an input cube.
- Parameters:
cube (numpy.ndarray) – 3D datacube (the final element contains the wavelength dimension)
wave (numpy.ndarray, optional) – 1D wavelength array. Only required if wavemin or wavemax are not None.
wavemin (float, optional) – Minimum wavelength (same units as wave) to be included in the whitelight image. You must provide wave as well if you want to reduce the wavelength range.
wavemax (float, optional) – Maximum wavelength (same units as wave) to be included in the whitelight image. You must provide wave as well if you want to reduce the wavelength range.
- Returns:
Whitelight image of the input cube.
- Return type:
- pypeit.core.datacube.make_whitelight_frompixels(all_ra, all_dec, all_wave, all_sci, all_wghts, all_idx, dspat, all_ivar=None, whitelightWCS=None, numra=None, numdec=None, trim=1)[source]
Generate a whitelight image using the individual pixels of every input frame
- Parameters:
all_ra (numpy.ndarray) – 1D flattened array containing the RA values of each pixel from all spec2d files
all_dec (numpy.ndarray) – 1D flattened array containing the DEC values of each pixel from all spec2d files
all_wave (numpy.ndarray) – 1D flattened array containing the wavelength values of each pixel from all spec2d files
all_sci (numpy.ndarray) – 1D flattened array containing the counts of each pixel from all spec2d files
all_wghts (numpy.ndarray) – 1D flattened array containing the weights attributed to each pixel from all spec2d files
all_idx (numpy.ndarray) – 1D flattened array containing an integer identifier indicating which spec2d file each pixel originates from. For example, a 0 would indicate that a pixel originates from the first spec2d frame listed in the input file. a 1 would indicate that this pixel originates from the second spec2d file, and so forth.
dspat (float) – The size of each spaxel on the sky (in degrees)
all_ivar (numpy.ndarray, optional) – 1D flattened array containing of the inverse variance of each pixel from all spec2d files. If provided, inverse variance images will be calculated and returned for each white light image.
whitelightWCS (astropy.wcs.WCS, optional) – The WCS of a reference white light image. If supplied, you must also supply numra and numdec.
numra (int, optional) – Number of RA spaxels in the reference white light image
numdec (int, optional) – Number of DEC spaxels in the reference white light image
trim (int, optional) – Number of pixels to grow around a masked region
- Returns:
two 3D arrays will be returned, each of shape [N, M, numfiles], where N and M are the spatial dimensions of the combined white light images. The first array is a white light image, and the second array is the corresponding inverse variance image. If all_ivar is None, this will be an empty array.
- Return type:
- pypeit.core.datacube.subpixellate(output_wcs, all_ra, all_dec, all_wave, all_sci, all_ivar, all_wghts, all_spatpos, all_specpos, all_spatid, tilts, slits, astrom_trans, bins, all_idx=None, spec_subpixel=10, spat_subpixel=10, debug=False)[source]
Subpixellate the input data into a datacube. This algorithm splits each detector pixel into multiple subpixels, and then assigns each subpixel to a voxel. For example, if
spec_subpixel = spat_subpixel = 10
, then each detector pixel is divided into \(10^2=100\) subpixels. Alternatively, when spec_subpixel = spat_subpixel = 1, this corresponds to the nearest grid point (NGP) algorithm.Important Note: If spec_subpixel > 1 or spat_subpixel > 1, the errors will be correlated, and the covariance is not being tracked, so the errors will not be (quite) right. There is a tradeoff one has to make between sampling and better looking cubes, versus no sampling and better behaved errors.
- Parameters:
output_wcs (astropy.wcs.WCS) – Output world coordinate system.
all_ra (numpy.ndarray) – 1D flattened array containing the right ascension of each pixel (units = degrees)
all_dec (numpy.ndarray) – 1D flattened array containing the declination of each pixel (units = degrees)
all_wave (numpy.ndarray) – 1D flattened array containing the wavelength of each pixel (units = Angstroms)
all_sci (numpy.ndarray) – 1D flattened array containing the counts of each pixel from all spec2d files
all_ivar (numpy.ndarray) – 1D flattened array containing the inverse variance of each pixel from all spec2d files
all_wghts (numpy.ndarray) – 1D flattened array containing the weights of each pixel to be used in the combination
all_spatpos (numpy.ndarray) – 1D flattened array containing the detector pixel location in the spatial direction
all_specpos (numpy.ndarray) – 1D flattened array containing the detector pixel location in the spectral direction
all_spatid (numpy.ndarray) – 1D flattened array containing the spatid of each pixel
tilts (numpy.ndarray, list) – 2D wavelength tilts frame, or a list of tilt frames (see all_idx)
slits (
SlitTraceSet
, list) – Information stored about the slits, or a list of SlitTraceSet (see all_idx)astrom_trans (
AlignmentSplines
, list) – A Class containing the transformation between detector pixel coordinates and WCS pixel coordinates, or a list of Alignment Splines (see all_idx)bins (tuple) – A 3-tuple (x,y,z) containing the histogram bin edges in x,y spatial and z wavelength coordinates
all_idx (numpy.ndarray, optional) – If tilts, slits, and astrom_trans are lists, this should contain a 1D flattened array, of the same length as all_sci, containing the index the tilts, slits, and astrom_trans lists that corresponds to each pixel. Note that, in this case all of these lists need to be the same length.
spec_subpixel (
int
, optional) – What is the subpixellation factor in the spectral direction. Higher values give more reliable results, but note that the time required goes as (spec_subpixel * spat_subpixel
). The default value is 5, which divides each detector pixel into 5 subpixels in the spectral direction.spat_subpixel (
int
, optional) – What is the subpixellation factor in the spatial direction. Higher values give more reliable results, but note that the time required goes as (spec_subpixel * spat_subpixel
). The default value is 5, which divides each detector pixel into 5 subpixels in the spatial direction.debug (bool) – If True, a residuals cube will be output. If the datacube generation is correct, the distribution of pixels in the residual cube with no flux should have mean=0 and std=1.
- Returns:
Three or four numpy.ndarray objects containing (1) the datacube generated from the subpixellated inputs, (2) the corresponding variance cube, (3) the corresponding bad pixel mask cube, and (4) the residual cube. The latter is only returned if debug is True.
- Return type: