Source¶
- class mpdaf.sdetect.Source(header, lines=None, mag=None, z=None, spectra=None, images=None, cubes=None, tables=None, mask_invalid=True, filename=None, default_size=None)[source]¶
Bases:
object
This class contains a Source object.
- Parameters:
- header
astropy.io.fits.Header
FITS header instance.
- lines
astropy.table.Table
List of lines.
- mag
astropy.table.Table
List of magnitudes.
- z
astropy.table.Table
List of redshifts.
- spectra
dict
Spectra dictionary, keys give origin of spectra (
'tot'
for total spectrum, TBC). Values areSpectrum
objects.- images
dict
Images dictionary, keys give filter names (
'MUSE_WHITE'
for white image, TBC). Values areImage
objects.- cubes
dict
Dictionary containing small data cubes. Keys gives a description of the cube. Values are
Cube
objects.- tables
dict
Dictionary containing tables. Keys give a description of each table. Values are
astropy.table.Table
objects.- mask_invalidbool
If True (default), iterate on all columns of all tables to mask invalid values (Inf, NaN, and -9999).
- default_size
float
Default size for image extraction, in arcseconds.
- header
Attributes Summary
Default size image extraction, in arcseconds.
Methods Summary
add_FSF
(cube[, fieldmap])Add FSF keywords from the cube FSF keywords.
add_attr
(key, value[, desc, unit, fmt])Add a new attribute for the current Source object.
add_comment
(comment, author[, date])Add a user comment to the FITS header of the Source object.
add_cube
(cube, name[, size, lbda, ...])Extract a cube centered on the source center and append it to the cubes dictionary.
add_history
(text[, author, date])Add a history to the FITS header of the Source object.
add_image
(image, name[, size, minsize, ...])Extract an image centered on the source center from the input image and append it to the images dictionary.
add_line
(cols, values[, units, desc, fmt, match])Add a line to the lines table.
add_mag
(band, m, errm)Add a magnitude value to the mag table.
add_narrow_band_image_lbdaobs
(cube, tag, lbda)Create narrow-band image around an observed wavelength value.
add_narrow_band_images
(cube, z_desc[, eml, ...])Create narrow-band images from a redshift value and a catalog of lines.
add_seg_images
([tags, DIR, del_sex, ...])Run SExtractor on all images to create segmentation maps.
add_table
(tab, name[, columns, select_in, ...])Append an astropy table to the tables dictionary.
add_white_image
(cube[, size, unit_size])Compute the white images from the MUSE data cube and appends it to the images dictionary.
add_z
(desc, z[, errz])Add a redshift value to the z table.
crack_z
([eml, nlines, cols, z_desc, zguess])Estimate the best redshift matching the list of emission lines.
extract_spectra
(cube[, obj_mask, sky_mask, ...])Extract spectra from a data cube.
find_intersection_mask
(seg_tags[, inter_mask])Use the list of segmentation maps to compute the instersection mask.
find_sky_mask
(seg_tags[, sky_mask])Loop over all segmentation images and use the region where no object is detected in any segmentation map as our sky image.
find_union_mask
(seg_tags[, union_mask])Use the list of segmentation maps to compute the union mask.
from_data
(ID, ra, dec, origin[, proba, ...])Source constructor from a list of data.
from_file
(filename[, ext, mask_invalid])Source constructor from a FITS file.
get_FSF
()Return the FSF model if available in the FITS header.
info
()Print information.
masked_invalid
([tables])Mask where invalid values occur (NaNs or infs or -9999 or '').
remove_attr
(key)Remove an Source attribute from the FITS header of the Source object.
show_ima
(ax, name[, showcenter, cuts, cmap])Show image.
show_rgb
(ax, names[, showcenter, cuts])Show RGB composite image.
show_spec
(ax, name[, cuts, zero, sky, lines])Display a spectra.
sort_lines
([nlines_max])Sort lines by flux in descending order.
write
(filename[, overwrite])Write the source object in a FITS file.
Attributes Documentation
- default_size¶
Default size image extraction, in arcseconds.
If not set, the size from the white-light image (MUSE_WHITE) is used.
Methods Documentation
- add_FSF(cube, fieldmap=None)[source]¶
Add FSF keywords from the cube FSF keywords.
For mosaic, when multiple fields are mixed with different FSF, the mean FSF is computed using a fieldmap.
The cube header must contain the FSFMODE keyword, and currently only the MOFFAT1 mode is handled. If there are multiple fields (keyword NFIELDS>1), then a fieldmap is required. It can be found either in a FIELDMAP extension in the cube file, or provided with the
fieldmap
parameter. Then the keywordsFSFxxBET
,FSFxxFWA
, andFSFxxFWB
where xx is the field index (00 for one field), are used to get the Moffat parameters, and are written in the source header.
- add_attr(key, value, desc=None, unit=None, fmt=None)[source]¶
Add a new attribute for the current Source object. This attribute will be saved as a keyword in the primary FITS header. This method could also be used to update a simple Source attribute that is saved in the pyfits header.
Equivalent to
self.key = (value, comment)
.- Parameters:
- key
str
Attribute name
- valueint/float/str
Attribute value
- desc
str
Attribute description
- unit
astropy.units.Unit
Attribute units
- fmt
str
Attribute format (‘.2f’ for example)
- key
- add_comment(comment, author, date=None)[source]¶
Add a user comment to the FITS header of the Source object.
- Parameters:
- comment
str
Comment.
- author
str
Initials of the author.
- date
datetime.date
Date, by default the current local date is used.
- comment
- add_cube(cube, name, size=None, lbda=None, add_white=False, unit_size=Unit('arcsec'), unit_wave=Unit('Angstrom'))[source]¶
Extract a cube centered on the source center and append it to the cubes dictionary.
Extracted cube saved in
self.cubes[name]
.- Parameters:
- cube
Cube
Input cube MPDAF object.
- name
str
Name used to distinguish this cube
- size
float
The size to extract. It corresponds to the size along the delta axis and the image is square. If None, the size of the white image extension is taken if it exists.
- lbda(
float
,float
) orNone
If not None, tuple giving the wavelength range.
- add_whitebool
Add white image from the extracted cube.
- unit_size
astropy.units.Unit
Unit of the size value (arcseconds by default). If None, size is in pixels.
- unit_wave
astropy.units.Unit
Wavelengths unit (angstrom by default). If None, inputs are in pixels.
- cube
- add_history(text, author='', date=None)[source]¶
Add a history to the FITS header of the Source object.
- Parameters:
- text
str
History text.
- author
str
Initials of the author.
- date
datetime.date
Date, by default the current local date is used.
- text
- add_image(image, name, size=None, minsize=2.0, unit_size=Unit('arcsec'), rotate=False, order=1)[source]¶
Extract an image centered on the source center from the input image and append it to the images dictionary.
Extracted image saved in
self.images[name]
.- Parameters:
- image
Image
Input image MPDAF object.
- name
str
Name used to distinguish this image
- size
float
The size to extract. It corresponds to the size along the delta axis and the image is square. If None, the size of the white image extension is taken if it exists.
- unit_size
astropy.units.Unit
Unit of
size
andminsize
. Arcseconds by default (use None for size in pixels).- minsize
float
The minimum size of the output image.
- rotatebool
if True, the image is rotated to the same PA as the white-light image.
- order
int
The order of the prefilter that is applied by the affine transform function for the rotation.
- image
- add_line(cols, values, units=None, desc=None, fmt=None, match=None)[source]¶
Add a line to the lines table.
- Parameters:
- cols
list
ofstr
Names of the columns
- valueslist<int/float/str>
List of corresponding values
- unitslist<astropy.units>
Unit of each column
- desc
list
ofstr
Description of each column
- fmt
list
ofstr
Format of each column.
- match(
str
, float/int/str, bool) Tuple (key, value, False/True) that gives the key to match the added line with an existing line. eg (‘LINE’,’LYALPHA1216’, True) If the boolean is True, the line will be added even if it is not matched.
- cols
- add_narrow_band_image_lbdaobs(cube, tag, lbda, size=None, unit_size=Unit('arcsec'), width=8, is_sum=False, subtract_off=True, margin=10.0, fband=3.0, median_filter=0, method='mean')[source]¶
Create narrow-band image around an observed wavelength value.
- Parameters:
- cube
Cube
MUSE data cube.
- tag
str
key used to identify the new narrow-band image in the images dictionary.
- lbda
float
Observed wavelength value in angstrom.
- size
float
The total size to extract in arcseconds. It corresponds to the size along the delta axis and the image is square. If None, the size of the white image extension is taken if it exists.
- unit_size
astropy.units.Unit
unit of the size value (arcseconds by default) If None, size is in pixels
- width
float
Angstrom total width
- is_sumbool
if True the image is computed as the sum over the wavelength axis, otherwise this is the average. Deprecated, use “sum” as aggregation method.
- subtract_offbool
If True, subtracting off nearby data.
- margin
float
This off-band is offseted by margin wrt narrow-band limit (in angstrom).
- fband
float
The size of the off-band is fband*narrow-band width (in angstrom).
- method
str
Name of the Cube method used to aggregate the data. This method must accept the axis=0 parameter and return an image. Example: mean, sum, max.
- cube
- add_narrow_band_images(cube, z_desc, eml=None, size=None, unit_size=Unit('arcsec'), width=8, is_sum=False, subtract_off=True, margin=10.0, fband=3.0, median_filter=0, method='mean')[source]¶
Create narrow-band images from a redshift value and a catalog of lines.
Algorithm from Jarle Brinchmann (jarle@strw.leidenuniv.nl)
Narrow-band images are saved in
self.images['NB_']
.- Parameters:
- cube
Cube
MUSE data cube.
- z_desc
str
Redshift description. The redshift value corresponding to this description will be used.
- eml
dict
{float:str
} Full catalog of lines Dictionary: key is the wavelength value in Angstrom, value is the name of the line. if None, the following catalog is used:
eml = {1216 : 'LYALPHA', 1908: 'SUMCIII1907', 3727: 'SUMOII3726', 4863: 'HBETA' , 5007: 'OIII5007', 6564: 'HALPHA'}
- size
float
The total size to extract. It corresponds to the size along the delta axis and the image is square. If None, the size of the white image extension is taken if it exists.
- unit_size
astropy.units.Unit
unit of the size value (arcseconds by default) If None, size is in pixels
- width
float
Narrow-band width(in angstrom).
- is_sumbool
if True the image is computed as the sum over the wavelength axis, otherwise this is the average. Deprecated, use “sum” as aggregation method.
- subtract_offbool
If True, subtracting off nearby data. The method computes the subtracted flux by using the algorithm from Jarle Brinchmann (jarle@strw.leidenuniv.nl):
# if method = "mean" sub_flux = mean(flux[lbda1-margin-fband*(lbda2-lbda1)/2: lbda1-margin] + flux[lbda2+margin: lbda2+margin+fband*(lbda2-lbda1)/2]) # or if method = "sum": sub_flux = sum(flux[lbda1-margin-fband*(lbda2-lbda1)/2: lbda1-margin] + flux[lbda2+margin: lbda2+margin+fband*(lbda2-lbda1)/2]) /fband # or if median_filter > 0: sub_flux = median_filter in the wavelength axis of flux
- margin
float
This off-band is offseted by margin wrt narrow-band limit (in angstrom).
- fband
float
The size of the off-band is
fband x narrow-band width
(in angstrom).- median_filter
float
size of the median filter for background estimation (if set to 0, the classical off band images are used )
- method
str
Name of the Cube method used to aggregate the data. This method must accept the axis=0 parameter and return an image. Example: mean, sum, max.
- cube
- add_seg_images(tags=None, DIR=None, del_sex=True, save_seg_table=False, outdir='./', debug=False)[source]¶
Run SExtractor on all images to create segmentation maps.
SExtractor will use the
default.nnw
,default.param
,default.sex
and*.conv
files present in the current directory. If not present default parameter files are created or copied from the directory given in input (DIR).Algorithm from Jarle Brinchmann (jarle@strw.leidenuniv.nl)
- Parameters:
- tags
list
ofstr
List of tags of selected images
- DIR
str
Directory that contains the configuration files of sextractor
- del_sexbool
If False, configuration files of sextractor are not removed.
- save_seg_tablebool
If True, segmentation table are saved in the source tables dict
- outdir
str
Name of directory where temporary files are copied and SExtractor run if None, a temporary directory with a unique name is created and deleted at the end
- debugbool
if True, the output of SExtractor is logged as DEBUG and the created and temporary directory files are not deleted
- tags
- add_table(tab, name, columns=None, select_in=None, margin=0, ra=None, dec=None, col_dist=None, col_edgedist=None, digit=3)[source]¶
Append an astropy table to the tables dictionary.
- Parameters:
- tab
astropy.table.Table
ormpdaf.sdetect.Catalog
Input Table object.
- name
str
Name used to distinguish this table.
- columns
list
ofstr
List of column names to select.
- select_in
str
Name of the image (available in the source) to use for the WCS selection.
- margin
int
Margin from the edges (pixels) for the WCS selection.
- ra
str
Name of the RA column (degrees) for WCS selection and distance.
- dec
str
Name of the DEC column (degrees) for WCS selection and distance.
- col_dist
str
Name of the column with the distance to the source in arcsec. If None distance is not computed. Defaults to DIST.
- col_edgedist
str
Name of the column with the distance to the image edges. If None distance is not computed.
- digit
int
Number of digits to round distances, defaults to 3.
- tab
- add_white_image(cube, size=5, unit_size=Unit('arcsec'))[source]¶
Compute the white images from the MUSE data cube and appends it to the images dictionary.
White image saved in self.images[‘MUSE_WHITE’].
- Parameters:
- cube
Cube
MUSE data cube.
- size
float
The total size to extract in arcseconds. It corresponds to the size along the delta axis and the image is square. By default 5x5arcsec
- unit_size
astropy.units.Unit
unit of the size value (arcseconds by default) If None, size is in pixels
- cube
- crack_z(eml=None, nlines=inf, cols=('LBDA_OBS', 'FLUX'), z_desc='EMI', zguess=None)[source]¶
Estimate the best redshift matching the list of emission lines.
Algorithm from Johan Richard (johan.richard@univ-lyon1.fr).
This method saves the redshift values in
self.z
and lists the detected lines inself.lines
.self.info()
could be used to print the results.- Parameters:
- eml
dict
{float:str
} Full catalog of lines to test redshift Dictionary: key is the wavelength value in Angtsrom, value is the name of the line. if None, the following catalog is used:
emlines = { 1215.67: 'LYALPHA1216' , 1550.0: 'CIV1550' , 1908.0: 'CIII]1909' , 2326.0: 'CII2326' , 3726.032: '[OII]3726' , 3728.8149: '[OII]3729' , 3798.6001: 'HTHETA3799' , 3834.6599: 'HETA3835' , 3869.0: '[NEIII]3869' , 3888.7: 'HZETA3889' , 3967.0: '[NEIII]3967' , 4102.0: 'HDELTA4102' , 4340.0: 'HGAMMA4340' , 4861.3198: 'HBETA4861' , 4959.0: '[OIII]4959' , 5007.0: '[OIII]5007' , 6548.0: '[NII6548]' , 6562.7998: 'HALPHA6563' , 6583.0: '[NII]6583' , 6716.0: '[SII]6716' , 6731.0: '[SII]6731' }
- nlines
int
estimated the redshift if the number of emission lines is inferior to this value
- cols(
str
,str
) tuple (wavelength column name, flux column name) Two columns of self.lines that will be used to define the emission lines.
- z_desc
str
Estimated redshift will be saved in self.z table under these name.
- zguess
float
Guess redshift. Test if this redshift is a match and fills the detected lines
- eml
- extract_spectra(cube, obj_mask='MASK_UNION', sky_mask='MASK_SKY', tags_to_try=('MUSE_WHITE', 'NB_LYALPHA', 'NB_HALPHA', 'NB_SUMOII3726'), skysub=True, psf=None, beta=None, lbda=None, apertures=None, unit_wave=Unit('Angstrom'))[source]¶
Extract spectra from a data cube.
This method extracts several spectra from a data cube and from a list of narrow-band images (to define spectrum extraction apertures). First, it computes a subcube that has the same size along the spatial axis as the mask image given by
obj_mask
.Then, the no-weighting spectrum is computed as the sum of the subcube weighted by the mask of the object and saved in
self.spectra['MUSE_TOT']
.The weighted spectra are computed as the sum of the subcube weighted by the corresponding narrow-band image. They are saved in
self.spectra[nb_ima] for nb_ima in tags_to_try
.For the weighted spectra, with the psf and narrow-band images, the optimal extraction algorithm for CCD spectroscopy Horne, K. 1986 is used. See
mpdaf.sdetect.compute_optimal_spectrum
for more detail.- If
psf
is True: The potential PSF weighted spectrum is computed as the sum of the subcube weighted by multiplication of the mask of the object and the PSF. It is saved in self.spectra[‘MUSE_PSF’]
- If
skysub
is True: The local sky spectrum is computed as the average of the subcube weighted by the sky mask image. It is saved in
self.spectra['MUSE_SKY']
The other spectra are computed on the sky-subtracted subcube and they are saved in
self.spectra['*_SKYSUB']
.
- Parameters:
- cube
Cube
Input data cube.
- obj_mask
str
Name of the image that contains the mask of the object.
- sky_mask
str
Name of the sky mask image.
- tags_to_try
list
ofstr
List of narrow-band images.
- skysubbool
If True, a local sky subtraction is done.
- psf
numpy.ndarray
The PSF to use for PSF-weighted extraction. This can be a vector of length equal to the wavelength axis to give the FWHM of the Gaussian or Moffat PSF at each wavelength (in arcsec) or a cube with the PSF to use. No PSF-weighted extraction by default.
- beta
float
ornone
If not none, the PSF is a Moffat function with beta value, else it is a Gaussian.
- lbda(
float
,float
) ornone
If not none, tuple giving the wavelength range.
- unit_wave
astropy.units.Unit
Wavelengths unit (angstrom by default) If None, inputs are in pixels
- apertures
list
offloat
List of aperture radii (arcseconds) for which a spectrum is extracted.
- cube
- If
- find_intersection_mask(seg_tags, inter_mask='MASK_INTER')[source]¶
Use the list of segmentation maps to compute the instersection mask.
Algorithm from Jarle Brinchmann (jarle@strw.leidenuniv.nl):
1- Select on each segmentation map the object at the centre of the map. (the algo supposes that each objects have different labels) 2- compute the intersection of these selected objects
- find_sky_mask(seg_tags, sky_mask='MASK_SKY')[source]¶
Loop over all segmentation images and use the region where no object is detected in any segmentation map as our sky image.
Algorithm from Jarle Brinchmann (jarle@strw.leidenuniv.nl)
- find_union_mask(seg_tags, union_mask='MASK_UNION')[source]¶
Use the list of segmentation maps to compute the union mask.
Algorithm from Jarle Brinchmann (jarle@strw.leidenuniv.nl):
1- Select on each segmentation map the object at the centre of the map. (the algo supposes that each objects have different labels) 2- compute the union of these selected objects
- classmethod from_data(ID, ra, dec, origin, proba=None, confid=None, extras=None, **kwargs)[source]¶
Source constructor from a list of data.
Additional parameters are passed to the
Source
constructor.- Parameters:
- ID
int
ID of the source
- ra
double
Right ascension in degrees
- dec
double
Declination in degrees
- origin
tuple
(str
,str
,str
,str
) 1- Name of the detector software which creates this object 2- Version of the detector software which creates this object 3- Name of the FITS data cube from which this object has been extracted. 4- Version of the FITS data cube
- proba
float
Detection probability
- confid
int
Expert confidence index
- extras
dict
{key:value
} ordict
{key: (value
,comment
)} Extra header keywords
- ID
- classmethod from_file(filename, ext=None, mask_invalid=True)[source]¶
Source constructor from a FITS file.
- show_ima(ax, name, showcenter=None, cuts=None, cmap='gray_r', **kwargs)[source]¶
Show image.
- Parameters:
- ax
matplotlib.axes._subplots.AxesSubplot
Matplotlib axis instance (eg ax = fig.add_subplot(2,3,1)).
- name
str
Name of image to display.
- showcenter(
float
,str
) radius in arcsec and color used to plot a circle around the center of the source.
- cuts(
float
,float
) Minimum and maximum values to use for the scaling.
- cmap
matplotlib.cm
Color map.
- kwargs
matplotlib.artist.Artist
kwargs can be used to set additional plotting properties.
- ax
- show_rgb(ax, names, showcenter=None, cuts=None, **kwargs)[source]¶
Show RGB composite image.
- Parameters:
- ax
matplotlib.axes._subplots.AxesSubplot
Matplotlib axis instance (eg ax = fig.add_subplot(2,3,1)).
- names[
str
,str
,str
] List of images corresponding to the blue, green and red filters.
- showcenter(
float
,str
) radius in arcsec and color used to plot a circle around the center of the source.
- cuts[(
float
,float
), (float
,float
), (float
,float
)] Minimum and maximum values to use for the scaling corresponding to the blue, green and red filters.
- kwargs
matplotlib.artist.Artist
kwargs can be used to set additional plotting properties.
- ax
- Returns:
- ax
matplotlib
AxesImage
Matplotlib axis instance.
- images_aligned
list
ofImage
The input images, but all aligned to that with the highest resolution.
- ax
- show_spec(ax, name, cuts=None, zero=False, sky=None, lines=None, **kwargs)[source]¶
Display a spectra.
- Parameters:
- ax
matplotlib.axes._subplots.AxesSubplot
Matplotlib axis instance (eg ax = fig.add_subplot(2,3,1)).
- name
str
Name of spectra to display.
- cuts(
float
,float
) Minimum and maximum values to use for the scaling.
- zero
float
If True, the 0 flux line is plotted in black.
- sky
Spectrum
Sky spectra to overplot (default None).
- lines
str
Name of a columns of the lines table containing wavelength values. If not None, overplot red vertical lines at the given wavelengths.
- kwargs
matplotlib.artist.Artist
kwargs can be used to set additional plotting properties.
- ax