ifscube package

Submodules

ifscube.channel_maps module

ifscube.channel_maps.channelmaps(cube, lambda0, vel_min, vel_max, channels=6, continuum_width=300, continuum_options=None, log_flux=False, angular_scale=None, scale_bar=None, north_arrow=None, lower_threshold=1e-16, plot_opts={}, fig_opts={}, width_space=None, height_space=None, text_color='black', stroke_color='white', color_bar=True, center_mark=True, screen=True, xy_coords=None)

Creates velocity channel maps from a data cube.

Parameters:
  • lambda0 (number) – Central wavelength of the desired spectral feature.
  • vel_min (number) – Mininum velocity in kilometers per second.
  • vel_max (number) – Maximum velocity in kilometers per second.
  • channels (int) – Number of channel maps to build.
  • continuum_width (number) – Width in wavelength units for the continuum evaluation window.
  • continuum_options (dict) – Dicitionary of options to be passed to the spectools.continuum function.
  • log_flux (bool) – If True, takes the base 10 logarithm of the fluxes.
  • lower_threshold (number) – Minimum emission flux for plotting, after subtraction of the continuum level. Spaxels with flux values below lowerThreshold will be masked in the channel maps.
  • angular_scale (number) – The angular pixel scale, in arcsec/pix. By default it is readen from the header keyword CD1_1.
  • scale_bar (dict) – Places a scale bar with the size ‘scale_size’ in the y,x position ‘scale_pos’, in the first panel, labeled with text ‘scale_tex’.
  • north_arrow (dict) – Places reference arrows where north PA is ‘north_pa’ and east is rotated 90 degrees counterclockwise (when ‘east_side’ is 1) or clockwise (when ‘east_side’ is -1). The arrows have origin at position ‘arrow_pos’.
  • color_bar (bool) – If True draws a colorbar.
  • center_mark (bool) – If True, evaluates the continuum centroid and marks it with ‘plus’ sign.
  • plot_opts (dict) – Dictionary of options to be passed to pcolormesh.
  • fig_opts (dict) – Options passed to pyplot.figure.
  • width_space (float) – Horizontal gap between channel maps.
  • height_space (float) – Vertical gap between channel maps.
  • text_color (matplotlib.color) – The color of the annotated texts specifying the velocity bin, the scale bar and scale text and the reference arrows and text.
  • stroke_color (matplotlib.color) – The color of the the thin stroke drawn around texts and lines to increase contrast when those symbols appear over image areas of similar color.
  • screen (bool) – If screen is True the channel maps are shown on screen.
  • xy_coords (tuple) – Tuple with (x, y) coordinates for the plots. If None coordinates will be automatically calculated based on the central pixel of the image.
ifscube.channel_maps.rgb_line_compose(cube, lambdas, velmin, velmax, channels=9, continuum_width=300, continuum_opts=None, logFlux=True, angScale=None, scaleBar={}, northArrow={}, lowerThreshold=1e-16, plot_opts={}, fig_opts={}, wspace=None, hspace=None, text_color='black', stroke_color='white', center_mark=True, bg_color=(0.25, 0.25, 0.25))

Creates velocity channel maps from a data cube for two or three different lines and compose them as an RGB image attributing each channel map to one of the RGB channels.

Parameters:
  • lambdas (list) – Central wavelength of the desired spectral feature.
  • vmin (number) – Mininum velocity in kilometers per second.
  • vmax (number) – Maximum velocity in kilometers per second.
  • channels (integer) – Number of channel maps to build.
  • continuum_width (number) – Width in wavelength units for the continuum evaluation window.
  • continuum_opts (dictionary) – Dicitionary of options to be passed to the spectools.continuum function.
  • logFlux (Boolean) – If True, takes the base 10 logarithm of the fluxes.
  • lowerThreshold (number) – Minimum emission flux for plotting, after subtraction of the continuum level. Spaxels with flux values below lowerThreshold will be masked in the channel maps.
  • angScale (number) – The angular pixel scale, in arcsec/pix. By default it is readen from the header keyword CD1_1.
  • scaleBar (dictionary) – Places a scale bar with the size ‘scale_size’ in the y,x position ‘scale_pos’, in the first panel, labeled with text ‘scale_tex’.
  • northArrow (dictionary) – Places reference arrows where north PA is ‘north_pa’ and east is rotated 90 degrees counterclockwise (when ‘east_side’ is 1) or clockwise (when ‘east_side’ is -1). The arrows have origin at position ‘arrow_pos’.
  • center_mark (Boolean) – If True, evaluates the continuum centroid and marks it with ‘plus’ sign.
  • plot_opts (dict) – Dictionary of options to be passed to pcolormesh.
  • fig_opts (dict) – Options passed to pyplot.figure.
  • wspace (number) – Horizontal gap between channel maps.
  • hspace (number) – Vertical gap between channel maps.
  • text_color (matplotlib color) – The color of the annotated text specifying the velocity bin.
  • stroke_color (matplotlibcolor) – The color of the the thin stroke drawn around texts and lines to increase contrast when those symbols appear over image areas of similar color.
  • bg_color (matplotlibcolor) – The color of the background area, where data is missing.

ifscube.cubetools module

Functions for the analysis of integral field spectroscopy.

Author: Daniel Ruschel Dutra Website: https://github.com/danielrd6/ifscube

class ifscube.cubetools.NanSolution

Bases: object

ppxf(sol, galaxy)
ifscube.cubetools.aperture_spectrum(arr, x0=None, y0=None, radius=3, combine='sum')
ifscube.cubetools.append_config(config_file: str, fit_file: str) → None

Appends the configuration parameters as FITS table to the fit output file.

Parameters:
  • config_file (str) – Name of the configuration file.
  • fit_file (str) – Name of the file onto which the configuration table will be appended.
Returns:

Return type:

None

ifscube.cubetools.bound_updater(p0, bound_range, bounds=None)
class ifscube.cubetools.gmosdc(*args, **kwargs)

Bases: object

ifscube.cubetools.nan_to_nearest(d)

Replaces nan values with the value of the nearest pixel.

Parameters:d (numpy.ndarray, numpy.ma.core.MaskedArray) – The array to have the nan values replaced. It d is a masked array, the masked values will also be replaced.
Returns:
  • g (numpy.ma.core.MaskedArray) – Output array.
  • Description
  • ———–
ifscube.cubetools.peak_spaxel(cube)
ifscube.cubetools.rebin(arr, xbin, ybin, combine='sum', mask=None)
ifscube.cubetools.resample_stellar(binned, original, normalization_wavelength=5500, span=30)

Resamples binned stellar spectra to each individual spaxel.

Parameters:
  • binned – Binned cube with stellar spectra.
  • original – Original data cube
  • normalization_wavelength – Wavelength coordinate to perform the normalization.
  • span – Width in wavelength of the normalization window.
Returns:

A data cube with the stellar spectra resampled for each spaxel.

Return type:

unbinned

ifscube.cubetools.scale_bounds(bounds, flux_sf)
ifscube.cubetools.wlprojection(arr, wl, wl0, fwhm=10, filtertype='box')

Writes a projection of the data cube along the wavelength coordinate, with the flux given by a given type of filter.

Parameters:
  • arr (np.ndarray) – Array to projected.
  • wl (np.ndarray) – Wavelength coordinates of the arr pixels.
  • wl0 (float) – Central wavelength at the rest frame.
  • fwhm (float) – Full width at half maximum. See ‘filtertype’.
  • filtertype (string) –

    Type of function to be multiplied by the spectrum to return the argument for the integral.

    ’box’ = Box function that is zero everywhere and 1
    between wl0-fwhm/2 and wl0+fwhm/2.
    ’gaussian’ = Normalized gaussian function with center at
    wl0 and sigma = fwhm/(2*sqrt(2*log(2)))
Returns:

outim – The integrated flux of the cube times the filter.

Return type:

numpy.ndarray

ifscube.cubetools_astropy module

ifscube.datacube module

class ifscube.datacube.Cube(fname: str = None, scidata: Union[str, int] = 'SCI', variance: Union[str, int] = None, flags: Union[str, int] = None, stellar: Union[str, int] = None, primary: Union[str, int] = 'PRIMARY', redshift: float = None, wcs_axis: int = 3, spatial_mask: Union[numpy.ndarray, str] = None, nan_spaxels: str = 'all')

Bases: ifscube.onedspec.Spectrum

A class for dealing with IFS data cubes.

aperture_spectrum(radius=1.0, x0=None, y0=None, flag_threshold=0.5, combine: str = 'sum')

Makes an aperture spectrum out of the data cube.

Parameters:
  • radius (float) – Radius of the virtual aperture in pixels.
  • y0 (x0,) – Central coordinates of the aperture in pixels. If both are set to None the center of the datacube will be used.
  • flag_threshold (float) – Amount of flagged pixels for the output spectrum to also be flagged in this pixel.
  • combine (str) – Combination type. Can be either ‘sum’, ‘mean’ or ‘median’.
Returns:

s – Spectrum object.

Return type:

ifscube.onedspec.Spectrum

channel_maps(*args, **kwargs)
continuum(write_fits=False, output_image=None, fitting_window=None, continuum_options=None)

Evaluates a polynomial continuum for the whole cube and stores it in self.cont.

Parameters:
  • write_fits (bool) – Write the output in a FITS file
  • output_image (str) – Name of the output FITS file.
  • fitting_window (list) – A list containing the starting and ending wavelengths.
  • continuum_options (dict) – Dictionary of continuum fitting options.
gaussian_smooth(sigma=2, write_fits=False, outfile=None)

Performs a spatial gaussian convolution on the data cube.

Parameters:
  • sigma (float) – Sigma of the gaussian kernel.
  • write_fits (bool) – Writes the output to a FITS file.
  • outfile (str) – Name of the output file.
Returns:

  • gdata (numpy.ndarray) – Gaussian smoothed data.
  • gvar (numpy.ndarray) – Smoothed variance.

See also

scipy.ndimage.gaussian_filter()

linefit(p0, write_fits=False, out_image=None, overwrite=False, individual_spec=None, refit=False, suffix=None, update_bounds=False, bound_range=0.1, spiral_loop=False, spiral_center=None, refit_radius=3.0, sig_threshold=0.0, par_threshold=0, verbose=False, debug: bool = False, **kwargs)

Fits a spectral feature with a gaussian function and returns a map of measured properties. This is a wrapper for the scipy minimize function that basically iterates over the cube, has a formula for the reduced chi squared, and applies an internal scale factor to the flux.

Parameters:
  • p0 (iterable) – Initial guess for the fitting funcion, consisting of a list of 3N parameters for N components of function. In the case of a gaussian fucntion, these parameters must be given as [amplitude0, center0, sigma0, amplitude1, center1, …].
  • write_fits (bool) – Writes the results in a FITS file.
  • out_image (string) – Name of the FITS file in which to write the results.
  • overwrite (bool) – Overwrites previously generated output file.
  • individual_spec (tuple, string or None) –

    Selects an individual spectrum from the data cube to be fit. It can be a tuple representing the horizontal and vertical coordinates of the spaxel or a string. Possible string values are:

    • ’peak’: Spaxel with the highest integrated flux.
    • ’cofm’: Spaxel at the center of mass of the integrated flux.

    If None all the spectra in the data cube will be fit.

  • refit (boolean) – Use parameters from nearby successful fits as the initial guess for the next fit.
  • suffix (string) – String to be appended to the end of the input file name. Only used when out_image is None.
  • update_bounds (boolean) – If using refit, update the bounds for the next fit.
  • bound_range (float) – Fractional difference for updating the bounds when using refit.
  • spiral_loop (boolean) – Begins the fitting with the central spaxel and continues spiraling outwards.
  • spiral_center (tuple) – Central coordinates for the beginning of the spiral given as a tuple of two coordinates (x0, y0).
  • refit_radius (float) – Spaxels within this radius will be considered in the reffiting process.
  • sig_threshold (float) – Fits which return par_threshold below this number of times the local noise will be set to np.nan. If set to 0 this criteria is ignored.
  • par_threshold (integer) – Parameter which must be above the noise threshold to be considered a valid fit.
  • verbose (integer) – Verbosity level.
  • debug (bool) – If there is an exception when performing the fit, prints the spaxel coordinates.
  • **kwargs – Additional arguments passed to ifscube.onedspec.Spectrum.linefit.
Returns:

sol – A data cube with the solution for each spectrum occupying the respective position in the image, and each position in the first axis giving the different parameters of the fit.

Return type:

numpy.ndarray

See also

ifscube.onedspec.Spectrum.linefit(), scipy.optimize.curve_fit(), scipy.optimize.leastsq()

static lineflux(amplitude, sigma)

Calculates the flux in a line given the amplitude and sigma of the gaussian function that fits it.

loadfit(fname)

Loads the result of a previous fit, and put it in the appropriate variables for the plotfit function.

Parameters:fname (string) – Name of the FITS file generated by gmosdc.linefit.
Returns:
Return type:Nothing.
peak_coords(wl_center, wl_width, center_type='peak_cen', spatial_center=None, spatial_width=10)

Returns the coordinates of the centroid or the peak of a 2D flux distrubution.

Parameters:
  • wl_center (number) – The central wavelenght of the spectral region over which the cube is to be inegrated.
  • wl_width (number) – The wavelenght width of the spectral region over which the cube is to be integrated.
  • center_type (string) – Type of centering algorithm emploied. Options are: centroid, peak and peak_cen. Where, ‘peak_cen’ returns the centroid on a box ‘spatial_width’ wide, centered on the pixel corresponding to the peak value, ‘centroid’ returns position of the centroid of the values in the region, and, ‘peak’ returns the pixel position of the maximum value in the region.
  • spatial_center (tuple) – Central position of the spatial region where the center is calculated.
  • spatial_width (number) – Side size of the square spatial region where the center is calculated.
plotfit(x=None, y=None, show=True, axis=None, output='stdout')

Plots the spectrum and features just fitted.

Parameters:
  • x (int) – Horizontal coordinate of the desired spaxel.
  • y (int) – Vertical coordinate of the desired spaxel.
  • show (bool) – Shows the plot.
  • axis (matplotlib.pyplot.Axes, optional) – Instance of Axes. If None a new instance will be created.
  • output (str) –

    Selects whether to print the fit results. Available options are:

    • ’stdout’: prints to standard output.
    • ’return’: returns the output as a string.
Returns:

Return type:

Nothing.

plotspec(x, y, show_noise: bool = False, noise_smooth: float = 0.0, ax: matplotlib.axes._axes.Axes = None)

Plots the spectrum at coordinates x,y.

Parameters:
  • x,y (numbers or iterables) – If x and y are numbers plots the spectrum at the specific spaxel. If x and y are two element tuples plots the average between x[0],y[0] and x[1],y[1]
  • show_noise (bool) – Displays the noise spectrum as a filled area.
  • noise_smooth (float) – Sigma of the gaussian kernel for the noise smoothing.
  • ax (matplotlib.pyplot.Axes, optional) – Axes instance in which to plot the spectrum. If None a new instance will be created.
Returns:

Return type:

Nothing.

snr_eval(wl_range=(6050, 6200), continuum_options=None)

Measures the signal to noise ratio (SNR) for each spectrum in a data cube, returning an image of the SNR.

This method evaluates the SNR for each spectrum in a data cube by measuring the residuals of a polynomial continuum fit. The function CONTINUUM of the SPECTOOLS package is used to provide the continuum, with zero rejection iterations and a 3 order polynomial.

Parameters:
  • wl_range (array like) – An array like object containing two wavelength coordinates that define the SNR window at the rest frame.
  • continuum_options (dictionary) – Options for the continuum fitting function.
Returns:

snr – Image of the SNR for each spectrum.

Return type:

numpy.ndarray

spatial_mask
spatial_rebin(xbin, ybin, combine='mean')

Spatial undersampling of the datacube.

Parameters:
  • xbin (int) – Size of the bin in the horizontal direction.
  • ybin (int) – Size of the bin in the vertical direction.
  • combine (str) –
    Type of spectral combination.
    • ’mean’: The spectral flux is averaged over the spatial bin.
    • ’sum’: The spectral flux is summed over the spatial bin.
Returns:

Return type:

None.

voronoi_binning(target_snr=10.0, write_fits=False, outfile=None, overwrite=False, plot=False, flag_threshold=0.5, **kwargs)

Applies Voronoi binning to the data cube, using Cappellari’s Python implementation.

Parameters:
  • target_snr (float) – Desired signal to noise ratio of the binned pixels
  • write_fits (boolean) – Writes a FITS image with the output of the binning.
  • plot (bool) – Plots the binning results.
  • outfile (string) – Name of the output FITS file. If ‘None’ then the name of the original FITS file containing the data cube will be used as a root name, with ‘.bin’ appended to it.
  • overwrite (boolean) – Overwrites files with the same name given in ‘outfile’.
  • flag_threshold (float) – Bins with less than this fraction of unflagged pixels will be flagged.
  • **kwargs (dict) – Arguments passed to voronoi_2d_binning.
Returns:

Return type:

Nothing.

Notes

The output file contains two tables which outline the tesselation process. These are stored in the extensions ‘VOR’ and ‘VORPLUS’.

wlprojection(wl0, fwhm, filtertype='box', writefits=False, outimage='outimage.fits')

Writes a projection of the data cube along the wavelength coordinate.

Parameters:
  • wl0 (float) – Central wavelength at the rest frame.
  • fwhm (float) – Full width at half maximum. See ‘filtertype’.
  • filtertype (string) –

    Type of function to be multiplied by the spectrum to return the argument for the integral. Should be one of

    • ’box’: Box function that is zero everywhere and 1 between \(\lambda_0 \pm {\rm FWHM}/2\)
    • ’gaussian’: Normalized gaussian function with center at
      \(\lambda_0\) and \(\sigma = {\rm FWHM}/(2\sqrt{2\log(2)})\)
  • writefits (bool) – Writes the output to a FITS file.
  • outimage (string) – Name of the output image
Returns:

Return type:

Nothing.

write(file_name: str, overwrite=False)
write_binnedspec(doppler_correction=False)

Writes only one spectrum for each bin in a FITS file.

Parameters:doppler_correction (bool) – Apply Doppler correction.

ifscube.diagnostics module

class ifscube.diagnostics.bpt(ha, n2, hb, o3)

Bases: object

See Baldwin, Philips & Terlevich 198?,
Kauffmann et al. 2003
kauffmann2003()
plot(ax=None, fig=None, xlim=(-1.5, 0.5), ylim=(-1.2, 1.5), **kwargs)
class ifscube.diagnostics.whan_diagram(wha, flux_ha, flux_n2)

Bases: object

See Cid Fernandes, R. et al 2011 MNRAS 413, 1687.

plot(ax=None, fig=None, text_opts={}, xlim=None, ylim=None, **kwargs)

ifscube.fitter module

ifscube.fitter.clear_lock(lock_name)
ifscube.fitter.dofit(file_name, line_fit_args, overwrite, cube_type, loading, fit_type, config_file_name, plot=False, lock=False)
ifscube.fitter.main(fit_type)
ifscube.fitter.make_lock(file_name)
ifscube.fitter.spectrum_fit(data: ifscube.onedspec.Spectrum, **line_fit_args)

ifscube.gmos module

class ifscube.gmos.Cube(*args, **kwargs)

Bases: ifscube.datacube.Cube

A class for dealing with data cubes, originally written to work with GMOS IFU.

ifscube.imgtools module

ifscube.imgtools.gauss2d(x, y, p)

Returns a 2-dimensional gaussian function, following the equation

Parameters:
  • y (x,) – The coordinates, or array of coordinates, of the point of interest.
  • p (list) – Parameters of the gaussian function in the following order: [a, b, x0, y0, sx, sy]
Returns:

g – 2D Gaussian.

Return type:

numpy.ndarray

Notes

g(x, y) = a * exp(-(x-x0)**2/(2 * sx**2) - (y-y0)**2/(2 * sy**2))

ifscube.imgtools.image2slit(data, x0=None, y0=None, pa=0, slit_width=1)

Produces a virtual slit image out of a 2D image.

Parameters:
  • data (numpy.ndarray) – Input image.
  • y0 (x0,) – Center coordinates of the virtual slit.
  • pa (float) – Position angle in degrees, with North being zero, and increasing towards the east.
  • slit_width (float) – Slit width in pixels.
Returns:

  • mask (numpy.ndarray) – Masked used to generate the virtual slit.
  • slit_x, slit_y (numpy.ndarray) – Horizontal and vertical coordinates of the slit.
  • slit_r (numpy.ndarray) – Distance from the slit center.
  • slit_z (numpy.ndarray) – Data values at the virtual slit positions.

ifscube.imgtools.match_mosaic(img, method='2dgaussian')

Attempts to match different parts of an image mosaic based on the assumption of a smooth distribution of surface brightness.

Parameters:
  • img (numpy.ndarray) – Array containing the mosaic to be matched
  • method (string) –
    Definition of the method to be employed. Possible options are:
    2dgaussian : Fits a 2D gaussian to the image and returns a
    ratio between the fit and the mosaic. This is best used when the source is known to be spatially unresolved, like the broad line region of an AGN, or a star.
Returns:

  • sol (numpy.ndarray) – Array containing the ratio image that best corrects for discrepancies in the flux of the mosaic.
  • r (minimize.optimize.OptimizeResult) – Optimization result.

See also

scipy.optimize.minimize()

ifscube.imgtools.pixel_distance(x, y, x0=0, y0=0)
ifscube.imgtools.rebin(arr, xbin, ybin, combine='sum', mask=None)

Re-samples images or data cubes spatially.

Parameters:
  • arr (numpy.ndarray) – Input data to be re-sampled.
  • ybin (xbin,) – Re-sampling factor in the horizontal and vertical directions.
  • combine (str) –
    Type of combination. Can be one of
    • sum : add the data in each bin
    • mean : averages the data in each bin
  • mask (numpy.ndarray) – Array of the same shape of arr in which masked elements are evaluated as True.
Returns:

new – Resulting re-sampled image or data cube.

Return type:

numpy.ndarray

ifscube.manga module

class ifscube.manga.IntegratedSpectrum(*args, **kwargs)

Bases: ifscube.onedspec.Spectrum

class ifscube.manga.cube(*args, **kwargs)

Bases: ifscube.datacube.Cube

ifscube.models module

class ifscube.models.DiskRotation(amplitude=100.0, c_0=1.0, p=1.25, phi_0=0.7854, theta=0.7854, v_sys=0.0, x_0=0.0, y_0=0.0, **kwargs)

Bases: astropy.modeling.core.Fittable2DModel

Observed radial velocity according to Bertola et al 1991.

Notes

This model returns a two-dimensional velocity field of observed radial velocities based on a rotation curve equal to

v_c(r) = A * r / ( (r ** 2 + c_0 ** 2) ** (p / 2) )

This equation was taken from Bertola et al. 1991 (ApJ, 373, 369B).

amplitude = Parameter('amplitude', value=100.0)
c_0 = Parameter('c_0', value=1.0)
static evaluate(x, y, amplitude, c_0, p, phi_0, theta, v_sys, x_0, y_0)
Parameters:
  • amplitude (float) – Integrated flux of the profile.
  • c_0 (float) – Concentration index.
  • p (float) – Velocity curve exponent.
  • phi_0 (float) – Position angle of the line of nodes in radians.
  • theta (float) – Inclination of the disk in radians(theta = 0 for a face-on disk).
Returns:

v – Velocity at coordinates x, y.

Return type:

numpy.ndarray

n_inputs = 2
n_outputs = 1
p = Parameter('p', value=1.25)
param_names = ('amplitude', 'c_0', 'p', 'phi_0', 'theta', 'v_sys', 'x_0', 'y_0')
phi_0 = Parameter('phi_0', value=0.7854)
theta = Parameter('theta', value=0.7854)
v_sys = Parameter('v_sys', value=0.0)
x_0 = Parameter('x_0', value=0.0)
y_0 = Parameter('y_0', value=0.0)
class ifscube.models.GaussHermite1D(amplitude=1, mean=0, stddev=1, h3=0, h4=0, **kwargs)

Bases: astropy.modeling.core.Fittable1DModel

Gauss-Hermite Series up to the fourth order.

Parameters:
  • amplitude (float) – Integrated flux of the profile.
  • mean (float) – Coordinate of the profile center.
  • stddev (float) – Sigma
  • h3 (float) – Coefficient of the third order element.
  • h4 (float) – Coefficient of the fourth order element.
  • Description
  • -----------
  • from Riffel, R. A. 2010 and van der Marel & Franx 1993. (Taken) –
amplitude = Parameter('amplitude', value=1.0)
static evaluate(x, amplitude, mean, stddev, h3, h4)

Evaluate the model on some input variables.

h3 = Parameter('h3', value=0.0)
h4 = Parameter('h4', value=0.0)
inputs = ('x',)
mean = Parameter('mean', value=0.0)
outputs = ('y',)
param_names = ('amplitude', 'mean', 'stddev', 'h3', 'h4')
stddev = Parameter('stddev', value=1.0)

ifscube.onedspec module

class ifscube.onedspec.Spectrum(fname: str = None, scidata: Union[str, int] = 'SCI', variance: Union[str, int] = None, flags: Union[str, int] = None, stellar: Union[str, int] = None, primary: Union[str, int] = 'PRIMARY', redshift: float = None, wcs_axis: int = None)

Bases: object

dn4000()

Notes

Dn4000 index based on Balogh et al. 1999 (ApJ, 527, 54).

The index is defined as the ratio between continuum fluxes at 3850A-3950A and 4000A-4100A.

red / blue

static doppler_correction(z, flux, wl)
plot(over_plot=True)

ifscube.parser module

class ifscube.parser.ConstraintParser(expr, parameter_names)

Bases: object

evaluate(parameter, scale_factor: float = 1.0, method: str = 'slsqp')
static isnumber(s)
tolist()
class ifscube.parser.LineFitParser(*args, **kwargs)

Bases: object

get_vars()
parse_line(line)

ifscube.plots module

ifscube.plots.velocity_width(results: dict)

ifscube.rotation module

class ifscube.rotation.Config(file_name)

Bases: configparser.ConfigParser

class ifscube.rotation.Rotation(input_data=None, extension=0, plane=0, model=None)

Bases: object

Class for fitting rotation models to velocity fields.

fit_model(maxiter=100)

Fits a rotation model to the data.

plot_results(contrast=1.0, contours=True, symmetric=True)

Plots the fit results.

print_solution()
update_bounds(bounds)
update_fixed(fixed)
update_model(parameters)

Updates the model parameters.

Parameters:parameters (configparser.Configparser) – Dictionary of parameters.
ifscube.rotation.main()

Fits a rotation model to a 2D array of velocities.

ifscube.spectools module

Spectools provide a number of small routines to work with 1D spectra in the form of numpy arrays or FITS files, depending on the function.

If the spectrum is specified as an array, the default is to assume that arr[:,0] are the wavelength coordinates and arr[:,1] are the flux points.

class ifscube.spectools.Constraints

Bases: object

static redshift(wla, wlb, rest0, rest1)
static same(ha, hb)
static same_intrinsic_sigma(ha, hb, instr_quadratic_difference=0)

Constraint for maintaining the same velocity dispersion.

Parameters:
  • ha (int) – Index for the sigma of the first spectral feature.
  • hb (int) – Index for the sigma of the first spectral feature.
  • instr_quadratic_difference (float) – (instr_ha ** 2) - (instr_hb ** 2). In km/s.
Returns:

d – Constraint dictionary for the SLSQP minimization method.

Return type:

dict

Notes

Same intrinsic velocity dispersion. sigma0_ha ** 2 = sigma0_hb ** 2 (sigma_ha ** 2 - inst_ha ** 2) = (sigma_hb ** 2 - inst_hb ** 2) (sigma_ha ** 2 - sigma_hb ** 2) - (inst_ha ** 2 - inst_hb ** 2) = 0

static sigma(sa, sb, wla, wlb)
ifscube.spectools.blackbody(x, t, coordinate='wavelength')

Evaluates the blackbody spectrum for a given temperature.

Parameters:
  • x (numpy.array) – Wavelength or frequency coordinates (CGS).
  • t (number) – Blackbody temperature in Kelvin
  • coordinate (string) – Specify the coordinates of x as ‘wavelength’ or ‘frequency’.
Returns:

b(x, T) – Flux density in cgs units.

Return type:

numpy.array

ifscube.spectools.closest(arr, value)

Returns the index of the array element that is numerically closest to the given value. The returned variable is an integer. This function expects a one-dimensional array.

>>> closest(np.arange(5), 3.7)
4
ifscube.spectools.continuum(x, y, output='ratio', degree=6, n_iterate=5, lower_threshold=2, upper_threshold=3, verbose=False, weights=None) → Union[Iterable, Callable]

Builds a polynomial continuum from segments of a spectrum, given in the form of wl and flux arrays.

Parameters:
  • x (array-like) – Independent variable
  • y (array-like) – y = f(x)
  • output (string) –

    Specifies what will be returned by the function

    ’ratio’ = ratio between fitted continuum and the spectrum ‘difference’ = difference between fitted continuum and the spectrum ‘function’ = continuum function evaluated at x

  • degree (integer) – Degree of polynomial for the fit
  • n_iterate (integer) – Number of rejection iterations
  • lower_threshold (float) – Lower threshold for point rejection in units of standard deviation of the residuals
  • upper_threshold (float) – Upper threshold for point rejection in units of standard deviation of the residuals
  • verbose (boolean) – Prints information about the fitting
  • weights (array-like) – Weights for continuum fitting. Must be the shape of x and y.
Returns:

c

c[0]: numpy.ndarray

Input x coordinates

c[1]: numpy.ndarray

See parameter “output”.

Return type:

tuple

ifscube.spectools.eqw(wl, flux, limits, continuum_iterate=5)

Measure the equivalent width of a feature in arr, defined by lims:

ifscube.spectools.feature_mask(wavelength, feature_wavelength, sigma, width=20, catch_error=False)

Creates the fit optimization window by selecting the appropriate portions of the spectrum.

Parameters:
  • wavelength (numpy.ndarray) – Wavelength coordinates
  • feature_wavelength (numpy.ndarray) – Central wavelength of the spectral features to be fit.
  • sigma (numpy.ndarray) – Sigma of each of the spectral features.
  • width (float) – Number of sigmas to each side that will be considered in the fit.
  • catch_error (bool) – Interrupts the program when the optimization window falls outside the available wavelengths.
Returns:

mask – True for all wavelengths within the window of interest for the fit.

Return type:

numpy.ndarray

ifscube.spectools.fit_gauss(x, y, p0=None, fit_center=True, fit_background=True)

Returns a gaussian function that fits the data. This function requires a background subtraction, meaning that the data should fall to zero within the fitting limits.

Parameters:
  • y (x,) – Independent variable and flux value vectors.
  • p0 (iterable) – Initial guesses for m, mu, sigma and bg respectively m = maximum value mu = center sigma = FWHM/(2*sqrt(2*log(2))) bg = background level
  • fit_center (boolean) – Chooses whether to fit center or accept the value given in p0[1].
  • fit_background (boolean) – Chooses whether to fit a horizontal background level or accept the value given in p0[3].
Returns:

ans – ans[0] = Lambda function with the fitted parameters ans[1] = Fit parameters

Return type:

tuple

ifscube.spectools.flags_to_mask(wavelength: numpy.ndarray, flags: numpy.ndarray) → list

Converts an array of flags into a list of masked intervals.

Parameters:
  • wavelength (numpy.ndarray) – Wavelength coordinates of the flags. These values will be the ones used in the mask list.
  • flags (numpy.ndarray) – Flags vector.
Returns:

new_mask – A list of wavelength pairs defining masked regions.

Return type:

list

ifscube.spectools.flambda2fnu(wl, fl)

Converts a flambda to fnu.

Parameters:
  • wl (1d-array) – Wavelength in angstroms
  • fl (1d-array) – Flux in ergs/s/cm^2/A
Returns:

fnu – Flux in Jy (10^23 erg/s/cm^2/Hz)

Return type:

1d-array

ifscube.spectools.fnu2flambda(fnu, wavelength)

Converts between flux units

Parameters:
  • fnu (number) – Flux density in W/m^2/Hz
  • wavelength (number) – Wavelength in microns
Returns:

  • flambda (number) – Flux density in W/m^2/um
  • Where does it come from? – d nu c
  • ——– = - ———-
  • d lambda lambda^2

ifscube.spectools.full_width_half_maximum(x: numpy.ndarray, y: numpy.ndarray, bg: list) → float

Evaluates the full width half maximum of y in units of x.

Parameters:
  • x (numpy.array) –
  • y (numpy.array) –
  • bg (list) – Background sampling limits
Returns:

full_width – Full width half maximum

Return type:

number

ifscube.spectools.get_wl(image, dimension=0, hdrext=0, dataext=0, dwlkey='CD1_1', wl0key='CRVAL1', pix0key='CRPIX1')

Obtains the wavelength coordinates from the header keywords of the FITS file image. The default keywords are CD1_1 for the delta lambda, CRVAL for the value of the first pixel and CRPIX1 for the number of the first pixel. These keywords are the standard for GEMINI images.

The function is prepared to work with Multi-Extesion FITS (MEF) files.

Parameters:
  • image (string) – Name of the FITS file containing the spectrum
  • dimension (integer) – Dimension of the dispersion direction
  • hdrext (number) – Extension that contains the header
  • dataext (number) – Extension that contains the actual spectral data
  • dwlkey (string) – Header keyword for the interval between two data points
  • wl0key (string) – Header keyword for the first pixel value
  • pix0key (string) – Header keyword for the first pixel coordinate
Returns:

wl – Wavelength coordinates for each data point

Return type:

numpy.array

ifscube.spectools.joinspec(x1, y1, x2, y2)

Joins two spectra

Parameters:
  • x1 (array) – Wavelength coordinates of the first spectrum
  • y1 (array) – Flux coordinates of the first spectrum
  • x2 (array) – Wavelength coordinates of the second spectrum
  • y2 (array) – Flux coordinates of the second spectrum
Returns:

  • x (array) – Joined spectrum wavelength coordinates
  • f (array) – Flux coordinates

ifscube.spectools.mask(arr, maskfile)

Eliminates masked points from a spectrum.

Parameters:
  • arr (ndarray) – Input spectrum for masking. The array can have any number of columns, provided that the wavelength is in the first one.
  • maskfile (string) – Name of the ASCII file containing the mask definitions, with one region per line defined by a lower and an upper limit in this order.
Returns:

b – An exact copy of arr without the points that lie within the regions defined by maskfile.

Return type:

ndarray

ifscube.spectools.natural_sort(sequence)
ifscube.spectools.normspec(x, y, wl, span)

Normalizes a copy of a 1D spectrum given as arr_in, with the wavelength as arr_in[:,0] and the flux as arr_in[:,1]. The normalization value is the average flux between wl-span/2. and wl+span/2..

Parameters:
  • x (numpy.array) – Input wavelength coordinates
  • y (numpy.array) – Input flux coordinates
  • wl (number) – Central wavelength for normalization
  • span (number) – Width of normalization window
ifscube.spectools.read_weights(wavelength: numpy.ndarray, file_name: str) → numpy.ndarray

Reads a weight definition mask from an input ASCII file and returns a vector of weights.

Parameters:
  • wavelength (np.ndarray) – Input wavelength coordinates.
  • file_name (str) –

    File containing the weight mask definition, where each line should have three values in the following order:

    <lower wavelength> <upper wavelength> <weight>

    By default every point has a weight of one.

Returns:

w – Vector of weights.

Return type:

np.ndarray

ifscube.spectools.remove_bg(x: numpy.ndarray, y: numpy.ndarray, sampling_limits: list, order: int = 1) → tuple

Removes a background function from one dimensional data.

Parameters:
  • x (numpy.array) –
  • y (numpy.array) –
  • sampling_limits (iterable) – A list of the background sampling limits
  • order (integer) – Polyfit order
Returns:

  • x (numpy.ndarray)
  • y_new (numpy.array)

ifscube.spectools.sigma_lambda(sigma_vel, rest_wl)
ifscube.spectools.spectrophotometry(spec: Callable, transmission: Callable, limits: Iterable, verbose: bool = False, get_filter_center: bool = False)

Evaluates the integrated flux for a given filter over a spectrum.

Parameters:
  • spec (function) – A function that describes the spectrum with only the wavelength or frequency as parameter.
  • transmission (function) – The same as the above but for the filter.
  • limits (iterable) – Lower and upper limits for the integration.
  • verbose (bool) – Verbosity.
  • get_filter_center (bool) – Returns the pivot wavelength.
Returns:

  • photometry (float) – Photometric data point that results from the integration and scaling of the filter over the spectrum (CGS).
  • Notes
  • ——
  • This method employs the scipy.integrate.quad function
  • to perform all the integrations, with a limit of 1000
  • and epsrel of 1.e-3.

ifscube.spectools.velocity_width(wavelength: numpy.ndarray, model: numpy.ndarray, data: numpy.ndarray, width: float = 80.0, smooth: float = None, clip_negative_flux: bool = True, sigma_factor: float = 5.0)

Evaluates the W80 parameter of a given emission fature.

Parameters:
  • wavelength (array-like) – Wavelength vector.
  • model (np.ndarray) – Modeled spectral feature.
  • data (np.ndarray) – Observed spectrum.
  • width (float) – Percentile velocity width. For instance width=80 for the W_80 index.
  • smooth (float) – Smoothing sigma to apply after the cumulative sum.
  • clip_negative_flux (bool) – Sets negative flux values to zero.
  • sigma_factor (float) – Radius of integration, from the modeled feature centroid, in units of the distance between the centroid and the 16 and 84 percentiles.
Returns:

res – Dictionary of results.

Return type:

dict

Notes

W80 is the width in velocity space which encompasses 80% of the light emitted in a given spectral feature. It is widely used as a proxy for identifying outflows of ionized gas in active galaxies. For instance, see Zakamska+2014 MNRAS.

ifscube.stats module

ifscube.stats.line_flux_error(flux, fwhm, sampling_interval, amplitude, amplitude_error)

Evaluates the uncertainty in the flux of the fitted line. Based on Lenz & Ayres 1992 and Wesson 2015

Parameters:
  • flux (float) – Integrated flux in the line.
  • fwhm (float) – Full width at half maximum.
  • sampling_interval (float) – Sampling step in wavelength, usually in angstroms per pixel.
  • amplitude (float) – The peak value of the flux.
  • amplitude_error (float) – Uncertainty in the peak flux.
Returns:

flux_error – Uncertainty in the integrated flux of the line.

Return type:

number

Module contents

class ifscube.Cube(fname: str = None, scidata: Union[str, int] = 'SCI', variance: Union[str, int] = None, flags: Union[str, int] = None, stellar: Union[str, int] = None, primary: Union[str, int] = 'PRIMARY', redshift: float = None, wcs_axis: int = 3, spatial_mask: Union[numpy.ndarray, str] = None, nan_spaxels: str = 'all')

Bases: ifscube.onedspec.Spectrum

A class for dealing with IFS data cubes.

aperture_spectrum(radius=1.0, x0=None, y0=None, flag_threshold=0.5, combine: str = 'sum')

Makes an aperture spectrum out of the data cube.

Parameters:
  • radius (float) – Radius of the virtual aperture in pixels.
  • y0 (x0,) – Central coordinates of the aperture in pixels. If both are set to None the center of the datacube will be used.
  • flag_threshold (float) – Amount of flagged pixels for the output spectrum to also be flagged in this pixel.
  • combine (str) – Combination type. Can be either ‘sum’, ‘mean’ or ‘median’.
Returns:

s – Spectrum object.

Return type:

ifscube.onedspec.Spectrum

channel_maps(*args, **kwargs)
continuum(write_fits=False, output_image=None, fitting_window=None, continuum_options=None)

Evaluates a polynomial continuum for the whole cube and stores it in self.cont.

Parameters:
  • write_fits (bool) – Write the output in a FITS file
  • output_image (str) – Name of the output FITS file.
  • fitting_window (list) – A list containing the starting and ending wavelengths.
  • continuum_options (dict) – Dictionary of continuum fitting options.
gaussian_smooth(sigma=2, write_fits=False, outfile=None)

Performs a spatial gaussian convolution on the data cube.

Parameters:
  • sigma (float) – Sigma of the gaussian kernel.
  • write_fits (bool) – Writes the output to a FITS file.
  • outfile (str) – Name of the output file.
Returns:

  • gdata (numpy.ndarray) – Gaussian smoothed data.
  • gvar (numpy.ndarray) – Smoothed variance.

See also

scipy.ndimage.gaussian_filter()

linefit(p0, write_fits=False, out_image=None, overwrite=False, individual_spec=None, refit=False, suffix=None, update_bounds=False, bound_range=0.1, spiral_loop=False, spiral_center=None, refit_radius=3.0, sig_threshold=0.0, par_threshold=0, verbose=False, debug: bool = False, **kwargs)

Fits a spectral feature with a gaussian function and returns a map of measured properties. This is a wrapper for the scipy minimize function that basically iterates over the cube, has a formula for the reduced chi squared, and applies an internal scale factor to the flux.

Parameters:
  • p0 (iterable) – Initial guess for the fitting funcion, consisting of a list of 3N parameters for N components of function. In the case of a gaussian fucntion, these parameters must be given as [amplitude0, center0, sigma0, amplitude1, center1, …].
  • write_fits (bool) – Writes the results in a FITS file.
  • out_image (string) – Name of the FITS file in which to write the results.
  • overwrite (bool) – Overwrites previously generated output file.
  • individual_spec (tuple, string or None) –

    Selects an individual spectrum from the data cube to be fit. It can be a tuple representing the horizontal and vertical coordinates of the spaxel or a string. Possible string values are:

    • ’peak’: Spaxel with the highest integrated flux.
    • ’cofm’: Spaxel at the center of mass of the integrated flux.

    If None all the spectra in the data cube will be fit.

  • refit (boolean) – Use parameters from nearby successful fits as the initial guess for the next fit.
  • suffix (string) – String to be appended to the end of the input file name. Only used when out_image is None.
  • update_bounds (boolean) – If using refit, update the bounds for the next fit.
  • bound_range (float) – Fractional difference for updating the bounds when using refit.
  • spiral_loop (boolean) – Begins the fitting with the central spaxel and continues spiraling outwards.
  • spiral_center (tuple) – Central coordinates for the beginning of the spiral given as a tuple of two coordinates (x0, y0).
  • refit_radius (float) – Spaxels within this radius will be considered in the reffiting process.
  • sig_threshold (float) – Fits which return par_threshold below this number of times the local noise will be set to np.nan. If set to 0 this criteria is ignored.
  • par_threshold (integer) – Parameter which must be above the noise threshold to be considered a valid fit.
  • verbose (integer) – Verbosity level.
  • debug (bool) – If there is an exception when performing the fit, prints the spaxel coordinates.
  • **kwargs – Additional arguments passed to ifscube.onedspec.Spectrum.linefit.
Returns:

sol – A data cube with the solution for each spectrum occupying the respective position in the image, and each position in the first axis giving the different parameters of the fit.

Return type:

numpy.ndarray

See also

ifscube.onedspec.Spectrum.linefit(), scipy.optimize.curve_fit(), scipy.optimize.leastsq()

static lineflux(amplitude, sigma)

Calculates the flux in a line given the amplitude and sigma of the gaussian function that fits it.

loadfit(fname)

Loads the result of a previous fit, and put it in the appropriate variables for the plotfit function.

Parameters:fname (string) – Name of the FITS file generated by gmosdc.linefit.
Returns:
Return type:Nothing.
peak_coords(wl_center, wl_width, center_type='peak_cen', spatial_center=None, spatial_width=10)

Returns the coordinates of the centroid or the peak of a 2D flux distrubution.

Parameters:
  • wl_center (number) – The central wavelenght of the spectral region over which the cube is to be inegrated.
  • wl_width (number) – The wavelenght width of the spectral region over which the cube is to be integrated.
  • center_type (string) – Type of centering algorithm emploied. Options are: centroid, peak and peak_cen. Where, ‘peak_cen’ returns the centroid on a box ‘spatial_width’ wide, centered on the pixel corresponding to the peak value, ‘centroid’ returns position of the centroid of the values in the region, and, ‘peak’ returns the pixel position of the maximum value in the region.
  • spatial_center (tuple) – Central position of the spatial region where the center is calculated.
  • spatial_width (number) – Side size of the square spatial region where the center is calculated.
plotfit(x=None, y=None, show=True, axis=None, output='stdout')

Plots the spectrum and features just fitted.

Parameters:
  • x (int) – Horizontal coordinate of the desired spaxel.
  • y (int) – Vertical coordinate of the desired spaxel.
  • show (bool) – Shows the plot.
  • axis (matplotlib.pyplot.Axes, optional) – Instance of Axes. If None a new instance will be created.
  • output (str) –

    Selects whether to print the fit results. Available options are:

    • ’stdout’: prints to standard output.
    • ’return’: returns the output as a string.
Returns:

Return type:

Nothing.

plotspec(x, y, show_noise: bool = False, noise_smooth: float = 0.0, ax: matplotlib.axes._axes.Axes = None)

Plots the spectrum at coordinates x,y.

Parameters:
  • x,y (numbers or iterables) – If x and y are numbers plots the spectrum at the specific spaxel. If x and y are two element tuples plots the average between x[0],y[0] and x[1],y[1]
  • show_noise (bool) – Displays the noise spectrum as a filled area.
  • noise_smooth (float) – Sigma of the gaussian kernel for the noise smoothing.
  • ax (matplotlib.pyplot.Axes, optional) – Axes instance in which to plot the spectrum. If None a new instance will be created.
Returns:

Return type:

Nothing.

snr_eval(wl_range=(6050, 6200), continuum_options=None)

Measures the signal to noise ratio (SNR) for each spectrum in a data cube, returning an image of the SNR.

This method evaluates the SNR for each spectrum in a data cube by measuring the residuals of a polynomial continuum fit. The function CONTINUUM of the SPECTOOLS package is used to provide the continuum, with zero rejection iterations and a 3 order polynomial.

Parameters:
  • wl_range (array like) – An array like object containing two wavelength coordinates that define the SNR window at the rest frame.
  • continuum_options (dictionary) – Options for the continuum fitting function.
Returns:

snr – Image of the SNR for each spectrum.

Return type:

numpy.ndarray

spatial_mask
spatial_rebin(xbin, ybin, combine='mean')

Spatial undersampling of the datacube.

Parameters:
  • xbin (int) – Size of the bin in the horizontal direction.
  • ybin (int) – Size of the bin in the vertical direction.
  • combine (str) –
    Type of spectral combination.
    • ’mean’: The spectral flux is averaged over the spatial bin.
    • ’sum’: The spectral flux is summed over the spatial bin.
Returns:

Return type:

None.

voronoi_binning(target_snr=10.0, write_fits=False, outfile=None, overwrite=False, plot=False, flag_threshold=0.5, **kwargs)

Applies Voronoi binning to the data cube, using Cappellari’s Python implementation.

Parameters:
  • target_snr (float) – Desired signal to noise ratio of the binned pixels
  • write_fits (boolean) – Writes a FITS image with the output of the binning.
  • plot (bool) – Plots the binning results.
  • outfile (string) – Name of the output FITS file. If ‘None’ then the name of the original FITS file containing the data cube will be used as a root name, with ‘.bin’ appended to it.
  • overwrite (boolean) – Overwrites files with the same name given in ‘outfile’.
  • flag_threshold (float) – Bins with less than this fraction of unflagged pixels will be flagged.
  • **kwargs (dict) – Arguments passed to voronoi_2d_binning.
Returns:

Return type:

Nothing.

Notes

The output file contains two tables which outline the tesselation process. These are stored in the extensions ‘VOR’ and ‘VORPLUS’.

wlprojection(wl0, fwhm, filtertype='box', writefits=False, outimage='outimage.fits')

Writes a projection of the data cube along the wavelength coordinate.

Parameters:
  • wl0 (float) – Central wavelength at the rest frame.
  • fwhm (float) – Full width at half maximum. See ‘filtertype’.
  • filtertype (string) –

    Type of function to be multiplied by the spectrum to return the argument for the integral. Should be one of

    • ’box’: Box function that is zero everywhere and 1 between \(\lambda_0 \pm {\rm FWHM}/2\)
    • ’gaussian’: Normalized gaussian function with center at
      \(\lambda_0\) and \(\sigma = {\rm FWHM}/(2\sqrt{2\log(2)})\)
  • writefits (bool) – Writes the output to a FITS file.
  • outimage (string) – Name of the output image
Returns:

Return type:

Nothing.

write(file_name: str, overwrite=False)
write_binnedspec(doppler_correction=False)

Writes only one spectrum for each bin in a FITS file.

Parameters:doppler_correction (bool) – Apply Doppler correction.