pypeit.core.datacube module

Module containing routines used by 3D datacubes.

pypeit.core.datacube.align_user_offsets(ifu_ra, ifu_dec, ra_offset, dec_offset)[source]

Align the RA and DEC of all input frames, and then manually shift the cubes based on user-provided offsets. The offsets should be specified in arcseconds, and the ra_offset should include the cos(dec) factor.

Parameters:
  • ifu_ra (numpy.ndarray) – A list of RA values of the IFU (one value per frame)

  • ifu_dec (numpy.ndarray) – A list of Dec values of the IFU (one value per frame)

  • ra_offset (numpy.ndarray) – A list of RA offsets to be applied to the input pixel values (one value per frame). Note, the ra_offset MUST contain the cos(dec) factor. This is the number of degrees on the sky that represents the telescope offset.

  • dec_offset (numpy.ndarray) – A list of Dec offsets to be applied to the input pixel values (one value per frame). This is the number of degrees on the sky that represents the telescope offset.

Returns:

A tuple containing a new set of RA and Dec offsets for each frame. Both arrays are of type numpy.ndarray, and are in units of degrees.

pypeit.core.datacube.check_inputs(list_inputs)[source]

This function checks the inputs to several of the cube building routines, and makes sure they are all consistent. Often, this is to make check if all inputs are lists of the same length, or if all inputs are 2D numpy.ndarray. The goal of the routine is to return a consistent set of lists of the input.

Parameters:

list_inputs (list) – A list of inputs to check.

Returns:

list_inputs – A list of inputs that have been checked for consistency.

Return type:

list

pypeit.core.datacube.compute_weights(raImg, decImg, waveImg, sciImg, ivarImg, slitidImg, all_wcs, all_tilts, all_slits, all_align, all_dar, ra_offsets, dec_offsets, whitelight_img, dspat, dwv, ra_min=None, ra_max=None, dec_min=None, dec_max=None, wave_min=None, wave_max=None, sn_smooth_npix=None, weight_method='auto')[source]

Calculate wavelength dependent optimal weights. The weighting is currently based on a relative \((S/N)^2\) at each wavelength

Parameters:
  • raImg (numpy.ndarray, list) – A list of 2D array containing the RA of each pixel, with shape (nspec, nspat)

  • decImg (numpy.ndarray, list) – A list of 2D array containing the Dec of each pixel, with shape (nspec, nspat)

  • waveImg (numpy.ndarray, list) – A list of 2D array containing the wavelength of each pixel, with shape (nspec, nspat)

  • sciImg (numpy.ndarray, list) – A list of 2D array containing the science image of each pixel, with shape (nspec, nspat)

  • ivarImg (numpy.ndarray, list) – A list of 2D array containing the inverse variance image of each pixel, with shape (nspec, nspat)

  • slitidImg (numpy.ndarray, list) – A list of 2D array containing the slit ID of each pixel, with shape (nspec, nspat)

  • all_wcs (astropy.wcs.WCS, list) – A list of WCS objects, one for each frame.

  • all_tilts (numpy.ndarray, list) – 2D wavelength tilts frame, or a list of tilt frames

  • all_slits (SlitTraceSet, list) – Information stored about the slits, or a list of SlitTraceSet objects

  • all_align (AlignmentSplines, list) – A Class containing the transformation between detector pixel coordinates and WCS pixel coordinates, or a list of Alignment Splines.

  • all_dar (DARcorrection, list) – A Class containing the DAR correction information, or a list of DARcorrection classes. If a list, it must be the same length as astrom_trans.

  • ra_offsets (float, list) – RA offsets for each frame in units of degrees

  • dec_offsets (float, list) – Dec offsets for each frame in units of degrees

  • whitelight_img (numpy.ndarray) – A 2D array containing a white light 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.

  • weight_method (str, optional) –

    Weight method to be used in sn_weights(). Options are 'auto', 'constant', 'uniform', 'wave_dependent', 'relative', or 'ivar'. The default is 'auto'. Behavior is as follows:

    • 'auto': Use constant weights if rms_sn < 3.0, otherwise

      use wavelength dependent.

    • 'constant': Constant weights based on rms_sn**2

    • 'uniform': Uniform weighting.

    • 'wave_dependent': Wavelength dependent weights will be

      used irrespective of the rms_sn ratio. This option will not work well at low S/N ratio although it is useful for objects where only a small fraction of the spectral coverage has high S/N ratio (like high-z quasars).

    • 'relative': Calculate weights by fitting to the ratio of

      spectra? Note, relative weighting will only work well when there is at least one spectrum with a reasonable S/N, and a continuum. RJC note - This argument may only be better when the object being used has a strong continuum + emission lines. The reference spectrum is assigned a value of 1 for all wavelengths, and the weights of all other spectra will be determined relative to the reference spectrum. This is particularly useful if you are dealing with highly variable spectra (e.g. emission lines) and require a precision better than ~1 per cent.

    • 'ivar': Use inverse variance weighting. This is not well

      tested and should probably be deprecated.

Returns:

all_wghts – Either a 2D numpy.ndarray or a list of 2D numpy.ndarray arrays containing the optimal weights of each pixel for all frames, with shape (nspec, nspat).

Return type:

list

pypeit.core.datacube.compute_weights_frompix(raImg, decImg, waveImg, sciImg, ivarImg, slitidImg, dspat, dwv, mnmx_wv, wghtsImg, all_wcs, all_tilts, all_slits, all_align, all_dar, ra_offsets, dec_offsets, ra_min=None, ra_max=None, dec_min=None, dec_max=None, wave_min=None, wave_max=None, sn_smooth_npix=None, weight_method='auto', reference_image=None, whitelight_range=None, specname='PYPSPEC')[source]

Calculate wavelength dependent optimal weights. The weighting is currently based on a relative \((S/N)^2\) at each wavelength. Note, this function first prepares a whitelight image, and then calls compute_weights() to determine the appropriate weights of each pixel.

Parameters:
  • raImg (numpy.ndarray, list) – A list of 2D array containing the RA of each pixel, with shape (nspec, nspat)

  • decImg (numpy.ndarray, list) – A list of 2D array containing the Dec of each pixel, with shape (nspec, nspat)

  • waveImg (numpy.ndarray, list) – A list of 2D array containing the wavelength of each pixel, with shape (nspec, nspat)

  • sciImg (numpy.ndarray, list) – A list of 2D array containing the science image of each pixel, with shape (nspec, nspat)

  • ivarImg (numpy.ndarray, list) – A list of 2D array containing the inverse variance image of each pixel, with shape (nspec, nspat)

  • slitidImg (numpy.ndarray, list) – A list of 2D array containing the slit ID of each pixel, with shape (nspec, nspat)

  • dspat (float) – The size of each spaxel on the sky (in degrees)

  • dwv (float) – The size of each wavelength pixel (in Angstroms)

  • mnmx_wv (numpy.ndarray) – The minimum and maximum wavelengths of every slit and frame. The shape is (Nframes, Nslits, 2), The minimum and maximum wavelengths are stored in the [:,:,0] and [:,:,1] indices, respectively.

  • wghtsImg (numpy.ndarray, list) – A list of 2D array containing the weights of each pixel, with shape (nspec, nspat)

  • all_wcs (astropy.wcs.WCS, list) – A list of WCS objects, one for each frame.

  • all_tilts (numpy.ndarray, list) – 2D wavelength tilts frame, or a list of tilt frames

  • all_slits (SlitTraceSet, list) – Information stored about the slits, or a list of SlitTraceSet objects

  • all_align (AlignmentSplines, list) – A Class containing the transformation between detector pixel coordinates and WCS pixel coordinates, or a list of Alignment Splines.

  • all_dar (DARcorrection, list) – A Class containing the DAR correction information, or a list of DARcorrection classes. If a list, it must be the same length as astrom_trans.

  • ra_offsets (float, list) – RA offsets for each frame in units of degrees

  • dec_offsets (float, list) – Dec offsets for each frame in units of degrees

  • ra_min (float, optional) – Minimum RA of the WCS (degrees)

  • ra_max (float, optional) – Maximum RA of the WCS (degrees)

  • dec_min (float, optional) – Minimum Dec of the WCS (degrees)

  • dec_max (float, optional) – Maximum Dec of the WCS (degrees)

  • wave_min (float, optional) – Minimum wavelength of the WCS (degrees)

  • wave_max (float, optional) – Maximum wavelength of the WCS (degrees)

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

  • weight_method (str, optional) –

    Weight method to be used in sn_weights(). Options are 'auto', 'constant', 'uniform', 'wave_dependent', 'relative', or 'ivar'. The default is 'auto'. Behavior is as follows:

    • 'auto': Use constant weights if rms_sn < 3.0, otherwise

      use wavelength dependent.

    • 'constant': Constant weights based on rms_sn**2

    • 'uniform': Uniform weighting.

    • 'wave_dependent': Wavelength dependent weights will be

      used irrespective of the rms_sn ratio. This option will not work well at low S/N ratio although it is useful for objects where only a small fraction of the spectral coverage has high S/N ratio (like high-z quasars).

    • 'relative': Calculate weights by fitting to the ratio of

      spectra? Note, relative weighting will only work well when there is at least one spectrum with a reasonable S/N, and a continuum. RJC note - This argument may only be better when the object being used has a strong continuum + emission lines. The reference spectrum is assigned a value of 1 for all wavelengths, and the weights of all other spectra will be determined relative to the reference spectrum. This is particularly useful if you are dealing with highly variable spectra (e.g. emission lines) and require a precision better than ~1 per cent.

    • 'ivar': Use inverse variance weighting. This is not well

      tested and should probably be deprecated.

  • reference_image (numpy.ndarray) – Reference image to use for the determination of the highest S/N spaxel in the image.

  • specname (str) – Name of the spectrograph

Returns:

weights – a 1D array the same size as all_sci, containing relative wavelength dependent weights of each input pixel.

Return type:

numpy.ndarray

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 (suffix ref). The grating correction is then evaluated at the wavelength array given by wave_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:

numpy.ndarray

pypeit.core.datacube.create_wcs(raImg, decImg, waveImg, slitid_img_gpm, dspat, dwave, ra_offsets=None, dec_offsets=None, ra_min=None, ra_max=None, dec_min=None, dec_max=None, wave_min=None, wave_max=None, reference=None, 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:
  • raImg ((numpy.ndarray, list):) – A list of 2D array containing the RA of each pixel, with shape (nspec, nspat)

  • decImg ((numpy.ndarray, list):) – A list of 2D array containing the Dec of each pixel, with shape (nspec, nspat)

  • (numpy.ndarray (waveImg) – A list of 2D array containing the wavelength of each pixel, with shape (nspec, nspat)

  • list) – A list of 2D array containing the wavelength of each pixel, with shape (nspec, nspat)

  • slitid_img_gpm ((numpy.ndarray, list):) – A list of 2D array containing the spat ID of each pixel, with shape (nspec, nspat). A value of 0 indicates that the pixel is not on a slit. All other values indicate the slit spatial ID.

  • dspat (float) – Spatial size of each square voxel (in arcsec). The default is to use the values in cubepar.

  • dwave (float) – Linear wavelength step of each voxel (in Angstroms)

  • ra_offsets (list, optional) – List of RA offsets for each frame (degrees)

  • dec_offsets (list, optional) – List of Dec offsets for each frame (degrees)

  • ra_min (float, optional) – Minimum RA of the WCS (degrees)

  • ra_max (float, optional) – Maximum RA of the WCS (degrees)

  • dec_min (float, optional) – Minimum Dec of the WCS (degrees)

  • dec_max (float, optional) – Maximum Dec of the WCS (degrees)

  • wave_min (float, optional) – Minimum wavelength of the WCS (degrees)

  • wave_max (float, optional) – Maximum wavelength of the WCS (degrees)

  • reference (str, optional) – Filename of a fits file that contains a WCS in the Primary HDU.

  • 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:

float

pypeit.core.datacube.extract_standard_spec(stdcube, subpixel=20)[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

Returns:

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 image

  • intflux (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:

numpy.ndarray

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:
  • crval (list) – 3 element list containing the [RA, DEC, WAVELENGTH] of the reference pixel

  • cdelt (list) – 3 element list containing the delta values of the [RA, DEC, WAVELENGTH]

  • equinox (float, optional) – Equinox of the WCS

Returns:

astropy WCS to be used for the combined cube

Return type:

astropy.wcs.WCS

pypeit.core.datacube.generate_cube_subpixel(output_wcs, bins, sciImg, ivarImg, waveImg, slitid_img_gpm, wghtImg, all_wcs, tilts, slits, astrom_trans, all_dar, ra_offset, dec_offset, spec_subpixel=5, spat_subpixel=5, slice_subpixel=5, overwrite=False, outfile=None, whitelight_range=None, debug=False)[source]

Save a datacube using the subpixel algorithm. Refer to the subpixellate() docstring for further details about this algorithm

Parameters:
  • output_wcs (astropy.wcs.WCS) – Output world coordinate system.

  • bins (tuple) – A 3-tuple (x,y,z) containing the histogram bin edges in x,y spatial and z wavelength coordinates

  • sciImg (numpy.ndarray, list) – A list of 2D array containing the counts of each pixel. If a list, the shape of each numpy array is (nspec, nspat).

  • ivarImg (numpy.ndarray, list) – A list of 2D array containing the inverse variance of each pixel. If a list, the shape of each numpy array is (nspec, nspat).

  • waveImg (numpy.ndarray, list) – A list of 2D array containing the wavelength of each pixel. If a list, the shape of each numpy array is (nspec, nspat).

  • slitid_img_gpm (numpy.ndarray, list) – A list of 2D array containing the slitmask of each pixel. If a list, the shape of each numpy array is (nspec, nspat). A zero value indicates that a pixel is either not on a slit or it is a bad pixel. All other values are the slit spatial ID number.

  • wghtImg (numpy.ndarray, list) – A list of 2D array containing the weights of each pixel to be used in the combination. If a list, the shape of each numpy array is (nspec, nspat).

  • all_wcs (astropy.wcs.WCS, list) – A list of astropy.wcs.WCS objects, one for each spec2d file.

  • tilts (list) – A list of numpy.ndarray objects, one for each spec2d file, containing the tilts of each pixel. The shape of each numpy array is (nspec, nspat).

  • slits (pypeit.slittrace.SlitTraceSet, list) – A list of pypeit.slittrace.SlitTraceSet objects, one for each spec2d file, containing the properties of the slit for each spec2d file

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

  • all_dar (DARcorrection, list) – A Class containing the DAR correction information, or a list of DARcorrection classes. If a list, it must be the same length as astrom_trans.

  • ra_offset (float, list) – A float or list of floats containing the RA offset of each spec2d file

  • dec_offset (float, list) – A float or list of floats containing the DEC offset of each spec2d file

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

  • slice_subpixel (int, optional) – What is the subpixellation factor in the slice direction. Higher values give more reliable results, but note that the time required goes as (slice_subpixel). The default value is 5, which divides each IFU slice into 5 subslices in the slice direction.

  • overwrite (bool, optional) – If True, the output cube will be overwritten.

  • outfile (str, optional) – Filename to be used to save 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.

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

Returns:

Four numpy.ndarray objects containing (1) the datacube generated from the subpixellated inputs. The shape of the datacube is (nwave, nspat1, nspat2). (2) the corresponding error cube (standard deviation). The shape of the error cube is (nwave, nspat1, nspat2). (3) the corresponding bad pixel mask cube. The shape of the bad pixel mask cube is (nwave, nspat1, nspat2). (4) a 1D array containing the wavelength at each spectral coordinate of the datacube. The shape of the wavelength array is (nwave,).

Return type:

tuple

pypeit.core.datacube.generate_image_subpixel(image_wcs, bins, sciImg, ivarImg, waveImg, slitid_img_gpm, wghtImg, all_wcs, tilts, slits, astrom_trans, all_dar, ra_offset, dec_offset, spec_subpixel=5, spat_subpixel=5, slice_subpixel=5, 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.

  • bins (tuple) – A 3-tuple (x,y,z) containing the histogram bin edges in x,y spatial and z wavelength coordinates

  • sciImg (numpy.ndarray, list) – A list of 2D science images, or a single 2D image containing the science data.

  • ivarImg (numpy.ndarray, list) – A list of 2D inverse variance images, or a single 2D image containing the inverse variance data.

  • waveImg (numpy.ndarray, list) – A list of 2D wavelength images, or a single 2D image containing the wavelength data.

  • slitid_img_gpm (numpy.ndarray, list) – A list of 2D slit ID images, or a single 2D image containing the slit ID data.

  • wghtImg (numpy.ndarray, list) – A list of 2D weight images, or a single 2D image containing the weight data.

  • all_wcs (astropy.wcs.WCS, list) – A list of WCS objects, or a single WCS object containing the WCS information of each image.

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

  • all_dar (DARcorrection, list) – A Class containing the DAR correction information, or a list of DARcorrection classes. If a list, it must be the same length as astrom_trans.

  • ra_offset (float, list) – The RA offset to apply to each image, or a list of RA offsets.

  • dec_offset (float, list) – The DEC offset to apply to each image, or a list of DEC offsets.

  • 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 * slice_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 * slice_subpixel). The default value is 5, which divides each detector pixel into 5 subpixels in the spatial direction.

  • slice_subpixel (int, optional) – What is the subpixellation factor in the slice direction. Higher values give more reliable results, but note that the time required goes as (spec_subpixel * spat_subpixel * slice_subpixel). The default value is 5, which divides each IFU slice into 5 subpixels in the slice 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:

numpy.ndarray

pypeit.core.datacube.get_output_filename(fil, par_outfile, combine, idx=1)[source]

Get the output filename of a datacube, given the input

Parameters:
  • fil (str) – The spec2d filename.

  • par_outfile (str) – The user-specified output filename (see cubepar[‘output_filename’])

  • combine (bool) – Should the input frames be combined into a single datacube?

  • idx (int, optional) – Index of filename to be saved. Required if combine=False.

Returns:

The output filename to use.

Return type:

str

pypeit.core.datacube.get_output_whitelight_filename(outfile)[source]

Given the output filename of a datacube, create an appropriate whitelight fits file name

Parameters:

outfile (str) – The output filename used for the datacube.

Returns:

A string containing the output filename to use for the whitelight image.

pypeit.core.datacube.get_whitelight_pixels(all_wave, all_slitid, min_wl, max_wl)[source]

Determine which pixels are included within the specified wavelength range

Parameters:
  • all_wave (numpy.ndarray, list) – List of numpy.ndarray wavelength images. The length of the list is the number of spec2d frames. Each element of the list contains a wavelength image that provides the wavelength at each pixel on the detector, with shape is (nspec, nspat).

  • all_slitid (numpy.ndarray, list) – List of numpy.ndarray slitid images. The length of the list is the number of spec2d frames. Each element of the list contains a slitid image that provides the slit number at each pixel on the detector, with shape (nspec, nspat).

  • min_wl (float) – Minimum wavelength to consider

  • max_wl (float) – Maximum wavelength to consider

Returns:

The first element of the tuple is a list of numpy.ndarray slitid images (or a single numpy.ndarray slitid image if only one spec2d frame is provided), shape is (nspec, nspat), where a zero value corresponds to an excluded pixel (either outside the desired wavelength range, a bad pixel, a pixel not on the slit). All other pixels have a value equal to the slit number. The second element of the tuple is the wavelength difference between the maximum and minimum wavelength in the desired wavelength range.

Return type:

tuple

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:

list

pypeit.core.datacube.load_imageWCS(filename, ext=0)[source]

Load an image and return the image and the associated WCS.

Parameters:
  • filename (str) – A fits filename of an image to be used when generating white light images. Note, the fits file must have a valid 3D WCS.

  • ext (bool, optional) – The extension that contains the image and WCS

Returns:

An numpy.ndarray with the 2D image data and a astropy.wcs.WCS with the image WCS.

Return type:

tuple

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:

numpy.ndarray

pypeit.core.datacube.make_sensfunc(ss_file, senspar, blaze_wave=None, blaze_spline=None, grating_corr=False)[source]

Generate the sensitivity function from a standard star DataCube.

Parameters:
  • ss_file (str) – The relative path and filename of the standard star datacube. It should be fits format, and for full functionality, should ideally of the form DataCube.

  • senspar (SensFuncPar) – The parameters required for the sensitivity function computation.

  • blaze_wave (numpy.ndarray, optional) – Wavelength array used to construct blaze_spline

  • blaze_spline (scipy.interpolate.interp1d, optional) – Spline representation of the reference blaze function (based on the illumflat).

  • grating_corr (bool, optional) – If a grating correction should be performed, set this variable to True.

Returns:

A mask of the good sky pixels (True = good)

Return type:

numpy.ndarray

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:

A whitelight image of the input cube (of type numpy.ndarray).

pypeit.core.datacube.set_voxel_sampling(spatscale, specscale, dspat=None, dwv=None)[source]

This function checks if the spatial and spectral scales of all frames are consistent. If the user has not specified either the spatial or spectral scales, they will be set here.

Parameters:
  • spatscale (numpy.ndarray) – 2D array, shape is (N, 2), listing the native spatial scales of N spec2d frames. spatscale[:,0] refers to the spatial pixel scale of each frame spatscale[:,1] refers to the slicer scale of each frame Each element of the array must be in degrees

  • specscale (numpy.ndarray) – 1D array listing the native spectral scales of multiple frames. The length of this array should be equal to the number of frames you are using. Each element of the array must be in Angstrom

  • dspat (float, optional) – Spatial scale to use as the voxel spatial sampling. If None, a new value will be derived based on the inputs

  • dwv (float, optional) – Spectral scale to use as the voxel spectral sampling. If None, a new value will be derived based on the inputs

Returns:

  • _dspat (float) – Spatial sampling

  • _dwv (float) – Wavelength sampling

pypeit.core.datacube.subpixellate(output_wcs, bins, sciImg, ivarImg, waveImg, slitid_img_gpm, wghtImg, all_wcs, tilts, slits, astrom_trans, all_dar, ra_offset, dec_offset, spec_subpixel=5, spat_subpixel=5, slice_subpixel=5, debug=False)[source]

Subpixellate the input data into a datacube. This algorithm splits each detector pixel into multiple subpixels and each IFU slice into multiple subslices. Then, the algorithm assigns each subdivided detector pixel to a voxel. For example, if spec_subpixel = spat_subpixel = slice_subpixel = 5, then each detector pixel is divided into \(5^3=125\) subpixels. Alternatively, when spec_subpixel = spat_subpixel = slice_subpixel = 1, this corresponds to the nearest grid point (NGP) algorithm.

Important Note: If spec_subpixel > 1 or spat_subpixel > 1 or slice_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.

  • bins (tuple) – A 3-tuple (x,y,z) containing the histogram bin edges in x,y spatial and z wavelength coordinates

  • sciImg (numpy.ndarray, list) – A list of 2D array containing the counts of each pixel. The shape of each 2D array is (nspec, nspat).

  • ivarImg (numpy.ndarray, list) – A list of 2D array containing the inverse variance of each pixel. The shape of each 2D array is (nspec, nspat).

  • waveImg (numpy.ndarray, list) – A list of 2D array containing the wavelength of each pixel. The shape of each 2D array is (nspec, nspat).

  • slitid_img_gpm (numpy.ndarray, list) – A list of 2D array containing the slitmask of each pixel. The shape of each 2D array is (nspec, nspat). A zero value indicates that a pixel is either not on a slit or it is a bad pixel. All other values are the slit spatial ID number.

  • wghtImg (numpy.ndarray, list) – A list of 2D array containing the weights of each pixel to be used in the combination. The shape of each 2D array is (nspec, nspat).

  • all_wcs (astropy.wcs.WCS, list) – A list of astropy.wcs.WCS objects, one for each spec2d file

  • tilts (list) – A list of numpy.ndarray objects, one for each spec2d file, containing the tilts of each pixel. The shape of each 2D array is (nspec, nspat).

  • slits (pypeit.slittrace.SlitTraceSet, list) – A list of pypeit.slittrace.SlitTraceSet objects, one for each spec2d file, containing the properties of the slit for each spec2d file

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

  • all_dar (DARcorrection, list) – A Class containing the DAR correction information, or a list of DARcorrection classes. If a list, it must be the same length as astrom_trans.

  • ra_offset (float, list) – A float or list of floats containing the RA offset of each spec2d file relative to the first spec2d file

  • dec_offset (float, list) – A float or list of floats containing the DEC offset of each spec2d file relative to the first spec2d file

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

  • slice_subpixel (int, optional) – What is the subpixellation factor in the slice direction. Higher values give more reliable results, but note that the time required goes as (slice_subpixel). The default value is 5, which divides each IFU slice into 5 subslices in the slice direction.

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

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:

tuple

pypeit.core.datacube.wcs_bounds(raImg, decImg, waveImg, slitid_img_gpm, ra_offsets=None, dec_offsets=None, ra_min=None, ra_max=None, dec_min=None, dec_max=None, wave_min=None, wave_max=None)[source]

Calculate the bounds of the WCS and the expected edges of the voxels, based on user-specified parameters or the extremities of the data. This is a convenience function that calls the core function in datacube.

Parameters:
  • raImg ((numpy.ndarray, list):) – A list of 2D array containing the RA of each pixel, with shape (nspec, nspat)

  • decImg ((numpy.ndarray, list):) – A list of 2D array containing the Dec of each pixel, with shape (nspec, nspat)

  • (numpy.ndarray (waveImg) – A list of 2D array containing the wavelength of each pixel, with shape (nspec, nspat)

  • list) – A list of 2D array containing the wavelength of each pixel, with shape (nspec, nspat)

  • slitid_img_gpm ((numpy.ndarray, list):) – A list of 2D array containing the spat ID of each pixel, with shape (nspec, nspat). A value of 0 indicates that the pixel is not on a slit. All other values indicate the slit spatial ID.

  • ra_offsets (list, optional) – A list of the RA offsets for each frame

  • dec_offsets (list, optional) – A list of the Dec offsets for each frame

  • ra_min (float, optional) – Minimum RA of the WCS

  • ra_max (float, optional) – Maximum RA of the WCS

  • dec_min (float, optional) – Minimum Dec of the WCS

  • dec_max (float, optional) – Maximum Dec of the WCS

  • wave_min (float, optional) – Minimum wavelength of the WCS

  • wave_max (float, optional) – Maximum wavelength of the WCS

Returns:

  • _ra_min (float) – Minimum RA of the WCS

  • _ra_max (float) – Maximum RA of the WCS

  • _dec_min (float) – Minimum Dec of the WCS

  • _dec_max (float) – Maximum Dec of the WCS

  • _wave_min (float) – Minimum wavelength of the WCS

  • _wave_max (float) – Maximum wavelength of the WCS