Cube

class interferopy.cube.Cube(filename: str)[source]

Bases: object

Cube object represents interferometric data such as a map (2D) or a cube (3D). Allows performing tasks such as aperture integration or spectrum extraction. It is instantiated by reading in a fits file.

Example

filename="cube.image.fits"
ra, dec, freq = (205.533741, 9.477317341, 222.54) # coord for the emission line
c=Cube(filename)
# spectrum from the circular aperture, with associated error
flux, err, _ = c.spectrum(ra=ra, dec=dec, radius=1.5, calc_error=True)
# curve of growth up to the maximum radius, in steps of one pixel, in a chosen frequency channel
radius, flux, err, _ = c.growing_aperture(ra=ra, dec=dec, freq=freq, maxradius=5, calc_error=True)

Attributes Summary

radesys

Reference frame.

rms

Rms (root mean square) per channel.

Methods Summary

aperture_r([ra, dec, freq, px, py, channel, ...])

Obtain integrated flux within a circular aperture as a function of radius.

aperture_value([ra, dec, freq, channel, ...])

Get an aperture integrated value at the given coord.

deltavel([reffreq])

Compute channel width in velocity units (km/s).

distance_grid(px, py)

Grid of distances from the chosen pixel (can be subpixel accuracy).

findclumps_1kernel(output_file[, ...])

FINDCLUMP(s) algorithm (Decarli+2014,Walter+2016).

findclumps_full(output_file[, kernels, ...])

Run the full findclumps search and analysis for different kernel sizes, on the positive and/or negative cube(s):

freq2pix([freq])

Get the channel number of requested frequency.

get_radesys()

Read celestial reference frame from header.

get_rms()

Calculate rms for each channel of the cube.

growing_aperture([ra, dec, freq, maxradius, ...])

Compute curve of growth at the given coordinate position in a circular aperture, growing up to the max radius.

im2d([ch, freq, function])

Get a 2D map.

im_mask_channels(channels_to_mask[, mask_value])

Mask specific channels_to_mask to mask_value.

im_mask_values([value_to_mask, mask_value])

Mask specific values_to_mask in the image to mask_value.

log(text)

Basic logger function to allow better functionality in the future development.

make_snr_cube([store, filename, overwrite])

Convert cube into a signal-to-noise cube (i.e., divide by rms), in-place.

pix2freq([pz])

Get the frequency of a given channel.

pix2radec([px, py])

Convert pixels coordinates into ra and dec.

profile_r([ra, dec, freq, px, py, channel, ...])

Obtain azimuthaly averaged profile as a function of radius.

radec2pix([ra, dec, integer])

Convert ra and dec coordinates into pixels to be used as im[px, py].

set_radesys(value)

set_rms(value)

single_pixel_value([ra, dec, freq, channel, ...])

Get a single pixel value at the given coord.

spectrum([ra, dec, radius, px, py, channel, ...])

Extract the spectrum (for 3D cube) or a single flux density value (for 2D map) at a given coord (ra, dec) integrated within a circular aperture of a given radius.

vels(reffreq)

Compute velocities of all cube channels for a given reference frequency.

write_fitsfile(filename[, overwrite])

Write the current head and im to a new fitsfile (useful if modified).

Attributes Documentation

radesys

Reference frame.

rms

Rms (root mean square) per channel.

Methods Documentation

aperture_r(ra: Optional[float] = None, dec: Optional[float] = None, freq: Optional[float] = None, px: Optional[int] = None, py: Optional[int] = None, channel: int = 0, maxradius: float = 1.0, binspacing: Optional[float] = None, bins: Optional[list] = None, calc_error: bool = False)[source]

Obtain integrated flux within a circular aperture as a function of radius. If freq is undefined, will return the spectrum of the 3D cube. Alias function. Check the “growing_aperture” method for details.

Units: ra[deg], dec[deg], maxradius[arcsec], binspacing[arcsec], freq[GHz]

Returns

(radius, flux) or (radius, flux, err) if calc_error

aperture_value(ra: Optional[float] = None, dec: Optional[float] = None, freq: Optional[float] = None, channel: Optional[int] = None, radius: float = 1.0, calc_error=False)[source]

Get an aperture integrated value at the given coord. Optionally, return associated error. Wrapper function for the “spectrum” method.

Parameters
  • ra – Right ascention in degrees. If None, assumes central pixel.

  • dec – Declination in degrees.

  • freq – Frequency in GHz. If None, computes whole spectrum.

  • channel – Channel index (alternative to freq).

  • radius – Aperture radius in arcsec.

  • calc_error – If true, returns also an error estimate (rms x sqrt(num of beams in aperture)).

Returns

value or (value, error)

deltavel(reffreq: Optional[float] = None)[source]

Compute channel width in velocity units (km/s).

Parameters

reffreq – Computed around specific velocity. If empty, use referent one from the header.

Returns

Channel width in km/s. Sign reflects how channels are ordered.

distance_grid(px: float, py: float) ndarray[source]

Grid of distances from the chosen pixel (can be subpixel accuracy). Uses small angle approximation (simple Pythagorean distances). Distances are measured in numbers of pixels.

Parameters
  • px – Index of x coord.

  • py – Index of y coord.

Returns

2D grid of distances in pixels, same shape as the cube slice

findclumps_1kernel(output_file, rms_region=0.25, minwidth=3, sextractor_param_file='', clean_tmp=True, negative=False, verbose=False)[source]

FINDCLUMP(s) algorithm (Decarli+2014,Walter+2016). Takes the cube image and outputs the 3d (x,y,wavelength) position of clumps of a minimum SN specified. Works by using a top-hat filter on a rebinned version of the datacube.

Parameters
  • output_file – relative/absolute path to the output catalogue

  • rms_region – Region to compute the rms noise [2x2 array in image pixel coord]. If None, takes the central 25% pixels (square)

  • minwidth – Number of channels to bin

  • sextractor_param_file – if ‘’, gets the default.sex file in interferopy

  • clean_tmp – Whether to remove or not the temporary files created by Sextractor

Returns

findclumps_full(output_file, kernels=array([3, 5, 7, 9, 11, 13, 15, 17, 19]), rms_region=0.25, sextractor_param_file='', clean_tmp=True, ncores=1, run_search=True, run_positive=True, run_negative=True, run_crop=True, SNR_min=3, delta_offset_arcsec=2, delta_freq=0.1, run_fidelity=True, fidelity_bins=array([0.0, 0.2, 0.4, 0.6, 0.8, 1.0, 1.2, 1.4, 1.6, 1.8, 2.0, 2.2, 2.4, 2.6, 2.8, 3.0, 3.2, 3.4, 3.6, 3.8, 4.0, 4.2, 4.4, 4.6, 4.8, 5.0, 5.2, 5.4, 5.6, 5.8, 6.0, 6.2, 6.4, 6.6, 6.8, 7.0, 7.2, 7.4, 7.6, 7.8, 8.0, 8.2, 8.4, 8.6, 8.8, 9.0, 9.2, 9.4, 9.6, 9.8, 10.0, 10.2, 10.4, 10.6, 10.8, 11.0, 11.2, 11.4, 11.6, 11.8, 12.0, 12.2, 12.4, 12.6, 12.8, 13.0, 13.2, 13.4, 13.6, 13.8, 14.0, 14.2, 14.4, 14.6, 14.8]), min_SN_fit=4.0, fidelity_threshold=0.5, verbose=True)[source]

Run the full findclumps search and analysis for different kernel sizes, on the positive and/or negative cube(s):

  1. Search (see findclumps_1kernel())

  2. Crop doubles and trim candidates above min_SNR (see tools.run_line_stats_sex and tools_crop_doubles)

  3. Determine the fidelity and produce catalog of candidates (see tools.fidelity_analysis)

To skip any of the steps, use the appropriate flags. This is intended to be used to re-run parts of the search and/or analysis after changing the parameters (e.g., adding different kernels and/or modifying the cropping criteria).

Parameters
  • output_file – relative/absolute path to the output catalogue

  • kernels – list of odd kernel widths to use

  • rms_region – Region to compute the rms noise [2x2 array in image pixel coord]. If None, takes the central 25% pixels (square)

  • sextractor_param_file – path to sextractor param file

  • clean_tmp – cleanup intermediate sextractor files (disable for debugging)

  • ncores – number of cores for multiprocessing (done over kernels and pos/neg search)

  • run_search – if False will skip the search step

  • run_positive – if False will skip the positive search

  • run_negative – if False will skip the negative search

  • run_crop – if False will skip the crop step

  • SNR_min – Minimum SN of peaks to retain in the crop

  • delta_offset_arcsec – maximum offset to match detections in the cube [arcsec]

  • delta_freq – maximum frequency offset to match detections in the cube [GHz]

  • run_fidelity – if False will skip the fidelity analysis

  • fidelity_bins – SN bins for the fidelity analysis

  • min_SN_fit – minimum SN for fidelity fit

  • fidelity_threshold – Fidelity threshold above which to select candidates

  • verbose – increase overall verbosity

freq2pix(freq: Optional[float] = None)[source]

Get the channel number of requested frequency.

Parameters

freq – Frequency in GHz.

Returns

Channel index.

get_radesys()[source]

Read celestial reference frame from header. If not found, it will assume ICRS and raise a warning.

N.B. This is rarely used. Findclumps needs it to create ds9 region files.

get_rms()[source]

Calculate rms for each channel of the cube. Can take some time to compute on large cubes.

Returns

single rms value if 2D, array odf rms values for each channel if 3D cube

growing_aperture(ra: ~typing.Optional[float] = None, dec: ~typing.Optional[float] = None, freq: ~typing.Optional[float] = None, maxradius=1.0, binspacing: ~typing.Optional[float] = None, bins: ~typing.Optional[list] = None, px: ~typing.Optional[int] = None, py: ~typing.Optional[int] = None, channel: int = 0, profile=False, calc_error=False) -> (<class 'numpy.ndarray'>, <class 'numpy.ndarray'>, <class 'numpy.ndarray'>, <class 'numpy.ndarray'>)[source]

Compute curve of growth at the given coordinate position in a circular aperture, growing up to the max radius. Coordinates can be given in degrees (ra, dec) or pixels (px, py). If no coordinates are given, the center of the map is assumed. If no channel is provided, the first one is assumed.

Parameters
  • ra – Right ascention in degrees.

  • dec – Declination in degrees.

  • freq – Frequency in GHz, takes precedence over channel param.

  • maxradius – Max radius for aperture integration in arcsec.

  • binspacing – Resolution of the growth flux curve in arcsec, default is one pixel size.

  • bins – Custom bins for curve growth (1D np array).

  • px – Right ascention pixel coord (alternative to ra).

  • py – Declination pixel coord (alternative to dec).

  • channel – Index of the cube channel to take (alternative to freq).

  • profile – If True, compute azimuthally averaged profile, if False, compute cumulative aperture values

  • calc_error – Set to False to skip error calculations, if the rms computation is slow or not necessary.

Returns

radius, flux, err, npix - all 1D numpy arrays: aperture radius, cumulative flux within it, associated Poissionain error (based on number of beams inside the aprture and the map rms), number of pixels

im2d(ch: ~typing.Optional[int] = None, freq: ~typing.Optional[float] = None, function=<function sum>)[source]

Get a 2D map. Convenience function to avoid indexing notation.

Provided with a single ch (or freq, but not both) it will return that channel. Alternatively, ch or freq can also be a list with the [start, stop] of a slice, to which function is then applied to make a collapsed image (default: np.sum).

Parameters
  • ch – channel index, alternatively: a list with a slice [start, stop]` to which function is applied

  • freq – channel freq, alternatively: a list with a slice [start, stop]` to which function is applied

  • function – function to apply when ch or freq is a slice (default: np.sum), should be a (numpy) function that takes the argument axis=-1 and aggregates over it

Returns

2D numpy array.

im_mask_channels(channels_to_mask, mask_value=nan)[source]

Mask specific channels_to_mask to mask_value.

Parameters
  • channels_to_mask – list of channels to be masked (python indexing)

  • mask_value – the value of pixels to be masked (default: np.nan)

im_mask_values(value_to_mask=None, mask_value=nan)[source]

Mask specific values_to_mask in the image to mask_value.

For example, if regions with no data have values equal to zero or another fill value they can be set to np.nan. This is required for get_rms() to work correctly.

Parameters
  • value_to_mask – value to set to mask_value, if None the value of the first pixel is used `

  • mask_value – the value of pixels to be masked (default: np.nan)

log(text: str)[source]

Basic logger function to allow better functionality in the future development. All class functions print info through this wrapper. Could be extended to provide different levels of info, timestamps, or logging to a file.

make_snr_cube(store: bool = False, filename: Optional[str] = None, overwrite: bool = False)[source]

Convert cube into a signal-to-noise cube (i.e., divide by rms), in-place.

Parameters
  • store – save resulting cube under filename (default: False)

  • filenamefilename (default:None appends _snr to self.filename

  • overwrite – overwrite file if it exists

pix2freq(pz: Optional[int] = None)[source]

Get the frequency of a given channel. If no channel is given, the center channel is assumed.

Parameters

pz – Channel index.

Returns

Frequency in GHz.

pix2radec(px=None, py=None)[source]

Convert pixels coordinates into ra and dec. If no coords are given, the center of the map is assumed.

Parameters
  • px – X-axis index. Or list.

  • py – Y-axis index. Or list.

Returns

Coords ra, dec in degrees

profile_r(ra: Optional[float] = None, dec: Optional[float] = None, freq: Optional[float] = None, px: Optional[int] = None, py: Optional[int] = None, channel: int = 0, maxradius: float = 1.0, binspacing: Optional[float] = None, bins: Optional[list] = None, calc_error: bool = False)[source]

Obtain azimuthaly averaged profile as a function of radius. Alias function. Check the “growing_aperture” method for details.

Units: ra[deg], dec[deg], maxradius[arcsec], binspacing[arcsec], freq[GHz]

Returns

(radius, profile) or (radius, profile err) if calc_error

radec2pix(ra=None, dec=None, integer=True)[source]

Convert ra and dec coordinates into pixels to be used as im[px, py]. If no coords are given, the center of the map is assumed.

Parameters
  • ra – Right ascention in degrees. Or list.

  • dec – Declination in degrees. Or list.

  • integer – Round the coordinate to an integer value (to be used as indices).

Returns

Coords x and y in pixels (0 based index).

set_radesys(value)[source]
set_rms(value)[source]
single_pixel_value(ra: Optional[float] = None, dec: Optional[float] = None, freq: Optional[float] = None, channel: Optional[int] = None, calc_error: bool = False)[source]

Get a single pixel value at the given coord. Optionally, return associated error. Wrapper function for the “spectrum” method. If None, assumes central pixel.

Parameters
  • ra – Right ascention in degrees.

  • dec – Declination in degrees.

  • freq – Frequency in GHz. If None, computes whole spectrum.

  • channel – Channel index (alternative to freq).

  • calc_error – If true, returns also an error estimate (rms).

Returns

value or (value, error)

spectrum(ra: ~typing.Optional[float] = None, dec: ~typing.Optional[float] = None, radius=0.0, px: ~typing.Optional[int] = None, py: ~typing.Optional[int] = None, channel: ~typing.Optional[int] = None, freq: ~typing.Optional[float] = None, calc_error=False) -> (<class 'numpy.ndarray'>, <class 'numpy.ndarray'>, <class 'numpy.ndarray'>, <class 'numpy.ndarray'>)[source]

Extract the spectrum (for 3D cube) or a single flux density value (for 2D map) at a given coord (ra, dec) integrated within a circular aperture of a given radius. Coordinates can be given in degrees (ra, dec) or pixels (px, py). If no coordinates are given, the center of the map is assumed. Single plane can be chosen instead of full spectrum by defining freq or channel. If no radius is given, a single pixel value is extracted (usual units Jy/beam), otherwise aperture integrated spectrum is extracted (usual units of Jy).

Note: use the freqs field (or velocities method) to get the x-axis values.

Parameters
  • ra – Right ascention in degrees.

  • dec – Declination in degrees.

  • radius – Circular aperture radius in arcsec.

  • px – Right ascention pixel coord.

  • py – Declination pixel coord.

  • channel – Force extracton in a single channel of provided index (instead of the full cube).

  • freq – Frequency in GHz (alternative to channel).

  • calc_error – Set to False to skip error calculations, if the rms computation is slow or not necessary.

Returns

flux, err, npix: flux (spectrum), error estimate, and number of pixels in the aperture

vels(reffreq: float)[source]

Compute velocities of all cube channels for a given reference frequency.

Parameters

reffreq – Reference frequency in GHz. If empty, use referent one from the header.

Returns

Velocities in km/s.

write_fitsfile(filename: str, overwrite=False)[source]

Write the current head and im to a new fitsfile (useful if modified). This is still experimental and may fail on certain headers.

Parameters
  • filename – Path string to the output file.

  • overwrite – False by default.