PixTable

class mpdaf.drs.PixTable(filename, xpos=None, ypos=None, lbda=None, data=None, dq=None, stat=None, origin=None, weight=None, primary_header=None, save_as_ima=True, wcs=Unit("pix"), wave=Unit("Angstrom"), unit_data=Unit("ct"))[source]

Bases: object

PixTable class.

This class manages input/output for MUSE pixel table files. The FITS file is opened with memory mapping. Just the primary header and table dimensions are loaded. The methods get_xpos, get_ypos, get_lambda, get_data, get_dq, get_stat and get_origin must be used to get columns data.

Parameters:

filename : str

The FITS file name. None by default.

Attributes

filename (str) The FITS file name. None if any.
primary_header (astropy.io.fits.Header) The primary header.
nrows (int) Number of rows.
nifu (int) Number of merged IFUs that went into this pixel table.
wcs (astropy.units.Unit) Type of spatial coordinates of this pixel table (u.pix, u.deg or u.rad)
wave (astropy.units.Unit) Type of spectral coordinates of this pixel table
ima (bool) If True, pixtable is saved as multi-extension FITS image instead of FITS binary table.

Attributes Summary

fluxcal If True, this pixel table was flux-calibrated.
projection Return the projection type.
skysub If True, this pixel table was sky-subtracted.

Methods Summary

copy() Copy PixTable object in a new one and returns it.
divide_slice_median(skyref, pixmask)
extract([filename, sky, lbda, ifu, sl, ...]) Extracts a subset of a pixtable using the following criteria:
extract_from_mask(mask) Return a new pixtable extracted with the given mask.
get_column(name[, ksel]) Load a column and return it.
get_data([ksel, unit]) Load the data column and return it.
get_dq([ksel]) Load the dq column and return it.
get_exp() Load the exposure numbers and return it as a column.
get_keyword(key) Return the keyword value corresponding to key, adding the keyword prefix ('HIERARCH ESO DRS MUSE PIXTABLE').
get_lambda([ksel, unit]) Load the lambda column and return it.
get_origin([ksel]) Load the origin column and return it.
get_pos_sky([xpos, ypos]) Return the absolute position on the sky in degrees/pixel.
get_stat([ksel, unit]) Load the stat column and return it.
get_weight([ksel]) Load the weight column and return it.
get_xpos([ksel, unit]) Load the xpos column and return it.
get_ypos([ksel, unit]) Load the ypos column and return it.
info() Print information.
mask_column([maskfile]) Compute the mask column corresponding to a mask file.
origin2coords(origin) Converts the origin value and returns (ifu, slice, ypix, xpix).
origin2ifu(origin) Converts the origin value and returns the ifu number.
origin2slice(origin) Converts the origin value and returns the slice number.
origin2xoffset(origin[, ifu, sli]) Converts the origin value and returns the x coordinates offset.
origin2xpix(origin[, ifu, sli]) Converts the origin value and returns the x coordinates.
origin2ypix(origin) Converts the origin value and returns the y coordinates.
reconstruct_det_image([xstart, ystart, ...]) Reconstructs the image on the detector from the pixtable.
reconstruct_det_waveimage() Reconstructs an image of wavelength values on the detector from the pixtable.
select_exp(exp, col_exp) Return a mask corresponding to given exposure numbers.
select_ifus(ifus[, origin]) Return a mask corresponding to given ifus.
select_lambda(lbda[, unit]) Return a mask corresponding to the given wavelength range.
select_sky(sky) Return a mask corresponding to the given aperture on the sky
select_slices(slices[, origin]) Return a mask corresponding to given slices.
select_stacks(stacks[, origin]) Return a mask corresponding to given stacks.
select_xpix(xpix[, origin]) Return a mask corresponding to given detector pixels.
select_ypix(ypix[, origin]) Return a mask corresponding to given detector pixels.
selfcalibrate([pixmask, corr_clip, logfile]) Correct the background level of the slices.
set_column(name, data[, ksel]) Set a column (or a part of it).
set_data(data[, ksel, unit]) Set data column (or a part of it).
set_dq(dq[, ksel]) Set dq column (or a part of it).
set_keyword(key, val) Set the keyword value corresponding to key, adding the keyword prefix ('HIERARCH ESO DRS MUSE PIXTABLE').
set_lambda(lbda[, ksel, unit]) Set lambda column (or a part of it).
set_origin(origin[, ksel]) Set origin column (or a part of it).
set_stat(stat[, ksel, unit]) Set stat column (or a part of it).
set_weight(weight[, ksel]) Set weight column (or a part of it).
set_xpos(xpos[, ksel, unit]) Set xpos column (or a part of it).
set_ypos(ypos[, ksel, unit]) Set ypos column (or a part of it).
sky_ref([pixmask, dlbda, nmax, nclip, ...]) Compute the reference sky spectrum using sigma clipped median.
subtract_slice_median(skyref, pixmask)
write(filename[, save_as_ima]) Save the object in a FITS file.

Attributes Documentation

fluxcal

If True, this pixel table was flux-calibrated.

projection

Return the projection type.

  • ‘positioned’ for positioned pixtables
  • ‘projected’ for reduced pixtables
skysub

If True, this pixel table was sky-subtracted.

Methods Documentation

copy()[source]

Copy PixTable object in a new one and returns it.

divide_slice_median(skyref, pixmask)[source]
extract(filename=None, sky=None, lbda=None, ifu=None, sl=None, xpix=None, ypix=None, exp=None, stack=None, method='and')[source]

Extracts a subset of a pixtable using the following criteria:

  • aperture on the sky (center, size and shape)
  • wavelength range
  • IFU numbers
  • slice numbers
  • detector pixels
  • exposure numbers
  • stack numbers

The arguments can be either single value or a list of values to select multiple regions.

Parameters:

filename : str

The FITS filename used to save the resulted object.

sky : (float, float, float, char)

(y, x, size, shape) extract an aperture on the sky, defined by a center (y, x) in degrees/pixel, a shape (‘C’ for circular, ‘S’ for square) and size (radius or half side length) in arcsec/pixels.

lbda : (float, float)

(min, max) wavelength range in angstrom.

ifu : int or list

IFU number.

sl : int or list

Slice number on the CCD.

xpix : (int, int) or list

(min, max) pixel range along the X axis

ypix : (int, int) or list

(min, max) pixel range along the Y axis

exp : list of int

List of exposure numbers

stack : list of int

List of stack numbers

method : ‘and’ or ‘or’

Logical operation used to merge the criteria

Returns:

out : PixTable

extract_from_mask(mask)[source]

Return a new pixtable extracted with the given mask.

Parameters:

mask : numpy.ndarray

Mask (array of bool).

Returns:

out : PixTable

get_column(name, ksel=None)[source]

Load a column and return it.

Parameters:

name : str or attribute

Name of the column.

ksel : output of np.where

Elements depending on a condition.

Returns:

out : numpy.array

get_data(ksel=None, unit=None)[source]

Load the data column and return it.

Parameters:

ksel : output of np.where

Elements depending on a condition.

unit : astropy.units.Unit

Unit of the returned data.

Returns:

out : numpy.array

get_dq(ksel=None)[source]

Load the dq column and return it.

Parameters:

ksel : output of np.where

Elements depending on a condition.

Returns:

out : numpy.array

get_exp()[source]

Load the exposure numbers and return it as a column.

Returns:out : numpy.memmap
get_keyword(key)[source]

Return the keyword value corresponding to key, adding the keyword prefix ('HIERARCH ESO DRS MUSE PIXTABLE').

Parameters:

key : str

Keyword.

Returns:

out : keyword value

get_lambda(ksel=None, unit=None)[source]

Load the lambda column and return it.

Parameters:

ksel : output of np.where

Elements depending on a condition.

unit : astropy.units.Unit

Unit of the returned data.

Returns:

out : numpy.array

get_origin(ksel=None)[source]

Load the origin column and return it.

Parameters:

ksel : output of np.where

Elements depending on a condition.

Returns:

out : numpy.array

get_pos_sky(xpos=None, ypos=None)[source]

Return the absolute position on the sky in degrees/pixel.

Parameters:

xpos : numpy.array

xpos values

ypos : numpy.array

ypos values

Returns:

xpos_sky, ypos_sky : numpy.array, numpy.array

get_stat(ksel=None, unit=None)[source]

Load the stat column and return it.

Parameters:

ksel : output of np.where

Elements depending on a condition.

unit : astropy.units.Unit

Unit of the returned data.

Returns:

out : numpy.array

get_weight(ksel=None)[source]

Load the weight column and return it.

Parameters:

ksel : output of np.where

Elements depending on a condition.

Returns:

out : numpy.array

get_xpos(ksel=None, unit=None)[source]

Load the xpos column and return it.

Parameters:

ksel : output of np.where

Elements depending on a condition.

unit : astropy.units.Unit

Unit of the returned data.

Returns:

out : numpy.array

get_ypos(ksel=None, unit=None)[source]

Load the ypos column and return it.

Parameters:

ksel : output of np.where

Elements depending on a condition.

unit : astropy.units.Unit

Unit of the returned data.

Returns:

out : numpy.array

info()[source]

Print information.

mask_column(maskfile=None)[source]

Compute the mask column corresponding to a mask file.

Parameters:

maskfile : str

Path to a FITS image file with WCS information, used to mask out bright continuum objects present in the FoV. Values must be 0 for the background and >0 for objects.

Returns:

out : mpdaf.drs.PixTableMask

origin2coords(origin)[source]

Converts the origin value and returns (ifu, slice, ypix, xpix).

Parameters:

origin : int

Origin value.

Returns:

out : (int, int, float, float)

origin2ifu(origin)[source]

Converts the origin value and returns the ifu number.

Parameters:

origin : int

Origin value.

Returns:

out : int

origin2slice(origin)[source]

Converts the origin value and returns the slice number.

Parameters:

origin : int

Origin value.

Returns:

out : int

origin2xoffset(origin, ifu=None, sli=None)[source]

Converts the origin value and returns the x coordinates offset.

Parameters:

origin : int

Origin value.

Returns:

out : float

origin2xpix(origin, ifu=None, sli=None)[source]

Converts the origin value and returns the x coordinates.

Parameters:

origin : int

Origin value.

Returns:

out : float

origin2ypix(origin)[source]

Converts the origin value and returns the y coordinates.

Parameters:

origin : int

Origin value.

Returns:

out : float

reconstruct_det_image(xstart=None, ystart=None, xstop=None, ystop=None)[source]

Reconstructs the image on the detector from the pixtable.

The pixtable must concerns only one IFU, otherwise an exception is raised.

Returns:out : Image
reconstruct_det_waveimage()[source]

Reconstructs an image of wavelength values on the detector from the pixtable. The pixtable must concerns only one IFU, otherwise an exception is raised.

Returns:out : Image
select_exp(exp, col_exp)[source]

Return a mask corresponding to given exposure numbers.

Parameters:

exp : list of int

List of exposure numbers

Returns:

out : array of bool

mask

select_ifus(ifus, origin=None)[source]

Return a mask corresponding to given ifus.

Parameters:

ifu : int or list

IFU number.

Returns:

out : array of bool

mask

select_lambda(lbda, unit=Unit("Angstrom"))[source]

Return a mask corresponding to the given wavelength range.

Parameters:

lbda : (float, float)

(min, max) wavelength range in angstrom.

unit : astropy.units.Unit

Unit of the wavelengths in input.

Returns:

out : array of bool

mask

select_sky(sky)[source]

Return a mask corresponding to the given aperture on the sky (center, size and shape)

Parameters:

sky : (float, float, float, char)

(y, x, size, shape) extract an aperture on the sky, defined by a center (y, x) in degrees/pixel, a shape (‘C’ for circular, ‘S’ for square) and size (radius or half side length) in arcsec/pixels.

Returns:

out : array of bool

mask

select_slices(slices, origin=None)[source]

Return a mask corresponding to given slices.

Parameters:

slices : list of int

Slice number on the CCD.

Returns:

out : array of bool

mask

select_stacks(stacks, origin=None)[source]

Return a mask corresponding to given stacks.

Parameters:

stacks : list of int

Stacks numbers (1,2,3 or 4)

Returns

——-

out : array of bool

mask

select_xpix(xpix, origin=None)[source]

Return a mask corresponding to given detector pixels.

Parameters:

xpix : list

[(min, max)] pixel range along the X axis

Returns:

out : array of bool

mask

select_ypix(ypix, origin=None)[source]

Return a mask corresponding to given detector pixels.

Parameters:

ypix : list

[(min, max)] pixel range along the Y axis

Returns:

out : array of bool

mask

selfcalibrate(pixmask=None, corr_clip=15.0, logfile=None)[source]

Correct the background level of the slices.

This requires a Pixtable that is not sky-subtracted. It uses the mean sky level as a reference, and compute a multiplicative correction to apply to each slice to bring its background level to the reference one.

A mask of sources can be provided. A PixTableMask must be computed from a 2D mask file with mask_column. This mask will be merged with the ‘DQ’ column from the pixtable.

The method works on wavelength bins, defined in mpdaf.drs.pixtable.SKY_SEGMENTS. These bins have been chosen so that their edges do not fall on a sky line, with a 200A to 300A width. They can also be modified (by setting a new list to this global variable).

Then, the algorithm can be summarized as follow:

- foreach lambda bin:
    - foreach ifu:
        - foreach slice:
            - foreach pixel:
                - compute the mean flux
            - compute the median flux of the slice
        - compute the median flux of the ifu
        - foreach slice:
            - if too few points, use the mean ifu flux
    - compute the total mean flux

    - foreach ifu:
        - foreach slice:
            - compute the correction: total_flux / slice_flux
        - compute the mean and stddev of the corrections
        - for slices where |correction - mean| > corr_clip*std_dev:
            - use the mean ifu correction: total_flux / ifu_flux

- foreach ifu:
    - foreach slice:
        - foreach lambda bin:
            - reject spikes in the correction curves, using
              a comparison with the correction from the next and
              previous lambda bin.
Parameters:

pixmask : mpdaf.drs.PixTableMask

Column corresponding to a mask file (previously computed by mask_column).

corr_clip : float

Clipping threshold for slice corrections in one IFU.

logfile : str

Path to a file to which the log will be written.

Returns:

out : mpdaf.drs.PixTableAutoCalib

set_column(name, data, ksel=None)[source]

Set a column (or a part of it).

Parameters:

name : str or attribute

Name of the column.

data : numpy.array

data values

ksel : output of np.where

Elements depending on a condition.

set_data(data, ksel=None, unit=None)[source]

Set data column (or a part of it).

Parameters:

data : numpy.array

data values

ksel : output of np.where

Elements depending on a condition.

unit : astropy.units.Unit

unit of the data column in input.

set_dq(dq, ksel=None)[source]

Set dq column (or a part of it).

Parameters:

dq : numpy.array

dq values

ksel : output of np.where

Elements depending on a condition.

set_keyword(key, val)[source]

Set the keyword value corresponding to key, adding the keyword prefix ('HIERARCH ESO DRS MUSE PIXTABLE').

Parameters:

key : str

Keyword.

val : str or int or float

Value.

set_lambda(lbda, ksel=None, unit=None)[source]

Set lambda column (or a part of it).

Parameters:

lbda : numpy.array

lbda values

ksel : output of np.where

Elements depending on a condition.

unit : astropy.units.Unit

unit of the lambda column in input.

set_origin(origin, ksel=None)[source]

Set origin column (or a part of it).

Parameters:

origin : numpy.array

origin values

ksel : output of np.where

Elements depending on a condition.

set_stat(stat, ksel=None, unit=None)[source]

Set stat column (or a part of it).

Parameters:

stat : numpy.array

stat values

ksel : output of np.where

Elements depending on a condition.

unit : astropy.units.Unit

unit of the stat column in input.

set_weight(weight, ksel=None)[source]

Set weight column (or a part of it).

Parameters:

weight : numpy.array

weight values

ksel : output of np.where

Elements depending on a condition.

set_xpos(xpos, ksel=None, unit=None)[source]

Set xpos column (or a part of it).

Parameters:

xpos : numpy.array

xpos values

ksel : output of np.where

Elements depending on a condition.

unit : astropy.units.Unit

unit of the xpos column in input.

set_ypos(ypos, ksel=None, unit=None)[source]

Set ypos column (or a part of it).

Parameters:

ypos : numpy.array

ypos values

ksel : output of np.where

Elements depending on a condition.

unit : astropy.units.Unit

unit of the ypos column in input.

sky_ref(pixmask=None, dlbda=1.0, nmax=2, nclip=5.0, nstop=2, lmin=None, lmax=None)[source]

Compute the reference sky spectrum using sigma clipped median.

Algorithm from Kurt Soto (kurt.soto@phys.ethz.ch)

Parameters:

pixmask : mpdaf.drs.PixTableMask

Column corresponding to a mask file (previously computed by mask_column).

dlbda : double

Wavelength step in angstrom

nmax : int

Maximum number of clipping iterations

nclip : float or (float,float)

Number of sigma at which to clip. Single clipping parameter or lower/upper clipping parameters

nstop : int

If the number of not rejected pixels is less than this number, the clipping iterations stop.

Returns:

out : Spectrum

subtract_slice_median(skyref, pixmask)[source]
write(filename, save_as_ima=True)[source]

Save the object in a FITS file.

Parameters:

filename : str

The FITS filename.

save_as_ima : bool

If True, pixtable is saved as multi-extension FITS image instead of FITS binary table.