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

numpy.ndarray

numpy.floating

The spectral blaze function

blaze_wave

numpy.ndarray

numpy.floating

Wavelength array of the spectral blaze function

bpm

numpy.ndarray

numpy.uint8

Bad pixel mask of the datacube (0=good, 1=bad)

flux

numpy.ndarray

numpy.floating

Flux datacube in units of counts/s/Ang/arcsec^2or 10^-17 erg/s/cm^2/Ang/arcsec^2

fluxed

bool

Boolean indicating if the datacube is fluxed.

sensfunc

numpy.ndarray

numpy.floating

Sensitivity function 10^-17 erg/(counts/cm^2)

sig

numpy.ndarray

numpy.floating

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

astropy.io.fits.Header

filename

Filename to use when loading from file

Type

str

spect_meta

Parsed meta from the header

Type

dict

spectrograph

Build from PYP_SPEC

Type

pypeit.spectrographs.spectrograph.Spectrograph

_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

list

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^2or 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 a DataContainer 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 set otype=float and otype=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 have otype that are tuple, int, float, numpy.integer, numpy.floating, numpy.ndarray, or astropy.table.Table objects. E.g., datamodel values for otype cannot be dict.

classmethod from_file(ifile)[source]

Over-load pypeit.datamodel.DataContainer.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 pypeit.datamodel.DataContainer.to_file() to deal with the header

Parameters
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.calc_grating_corr(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

grat_corr (numpy.ndarray)

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 files

  • opts (dict) – coadd2d options associated with each spec2d file

  • spectrograph (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 (see default_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

numpy.ndarray

pypeit.core.datacube.dar_correction(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 location

  • temperature (float) – Outside ambient air temperature at location

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

chisq (float)

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

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

Return type

wave (numpy.ndarray)

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

gtwod (numpy.ndarray)

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 cube

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

gtwod (numpy.ndarray)

pypeit.core.datacube.generate_WCS(crval, cdelt, equinox=2000.0, name='Instrument Unknown')[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) – Equinox of the WCS

Returns

astropy WCS to be used for the combined cube

Return type

astropy.wcs.wcs.WCS

pypeit.core.datacube.generate_cube_ngp(outfile, hdr, all_sci, all_ivar, all_wghts, vox_coord, bins, overwrite=False, blaze_wave=None, blaze_spec=None, fluxcal=False, sensfunc=None, specname='PYP_SPEC', debug=False)[source]

Save a datacube using the Nearest Grid Point (NGP) algorithm.

Parameters
  • outfile (str) – Filename to be used to save the datacube

  • hdr (astropy.io.fits.header_) – Header of the output datacube (must contain WCS)

  • 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

  • vox_coord (numpy.ndarray) – The voxel coordinates of each pixel in the spec2d frames. vox_coord is returned by the function astropy.wcs.WCS.wcs_world2pix_ once a WCS is setup and every spec2d detector pixel has an RA, DEC, and WAVELENGTH.

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

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

  • blaze_wave (numpy.ndarray) – Wavelength array of the spectral blaze function

  • blaze_spec (numpy.ndarray) – Spectral blaze function

  • fluxcal (bool) – Are the data flux calibrated?

  • sensfunc (numpy.ndarray, None) – Sensitivity function that has been applied to the datacube

  • specname (str) – Name of the spectrograph

  • debug (bool) – Debug the code by writing out a residuals cube?

pypeit.core.datacube.generate_cube_subpixel(outfile, output_wcs, all_sci, all_ivar, all_wghts, all_wave, tilts, slits, slitid_img_gpm, astrom_trans, bins, spec_subpixel=10, spat_subpixel=10, overwrite=False, blaze_wave=None, blaze_spec=None, fluxcal=False, sensfunc=None, specname='PYP_SPEC', debug=False)[source]

Save a datacube using the subpixel algorithm. 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
  • outfile (str) – Filename to be used to save the datacube

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

  • 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_wave (numpy.ndarray) – 1D flattened array containing the wavelength of each pixel (units = Angstroms)

  • tilts (numpy.ndarray) – 2D wavelength tilts frame

  • slits (pypeit.slittrace.SlitTraceSet) – Information stored about the slits

  • slitid_img_gpm (numpy.ndarray) – An image indicating which pixels belong to a slit (0 = not on a slit or a masked pixel). Any positive value indicates the spatial ID of the pixel.

  • astrom_trans (pypeit.alignframe.AlignmentSplines) – A Class containing the transformation between detector pixel coordinates and WCS pixel coordinates

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

  • 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: erg/s/cm^2/Angstrom/arcsec^2 multiplied by the PYPEIT_FLUX_SCALE. Otherwise, the units are: counts/s/Angstrom/arcsec^2”)

  • sensfunc (numpy.ndarray, None, optional) – Sensitivity function that has been applied to the datacube

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

outfile (str)

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

gpm (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, None) – 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, None) – 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

wl_img (numpy.ndarray)

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

tuple

pypeit.core.datacube.make_whitelight_fromref(all_ra, all_dec, all_wave, all_sci, all_wghts, all_idx, dspat, ref_filename)[source]

Generate a whitelight image using the individual pixels of every input frame, based on a reference image. Note the, the reference image must have a well-defined WCS.

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)

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

Returns

two numpy.ndarray and one WCS will be returned. The first is a 2D reference image loaded from ref_filename. The second element is a 3D array of shape [N, M, numfiles], where N and M are the spatial dimensions of the combined white light images. The third is the WCS of the white light image.

Return type

tuple

pypeit.core.datacube.rebinND(img, shape)[source]

Rebin a 2D image to a smaller shape. For example, if img.shape=(100,100), then shape=(10,10) would take the mean of the first 10x10 pixels into a single output pixel, then the mean of the next 10x10 pixels will be output into the next pixel

Parameters
  • img (numpy.ndarray) – A 2D input image

  • shape (tuple) – The desired shape to be returned. The elements of img.shape should be an integer multiple of the elements of shape.

Returns

The input image rebinned to shape

Return type

img_out (numpy.ndarray)