Bispectrum

class turbustat.statistics.Bispectrum(img)[source] [edit on github]

Bases: turbustat.statistics.base_statistic.BaseStatisticMixIn

Computes the bispectrum (three-point correlation function) of the given image (Burkhart et al., 2010). The bispectrum and the bicoherence are returned. The bicoherence is a normalized version (real and to unity) of the bispectrum.

Parameters:
img : numpy.ndarray or astropy.io.fits.PrimaryHDU or astropy.io.fits.ImageHDU or spectral_cube.Projection or spectral_cube.Slice

2D image.

Examples

>>> from turbustat.statistics import Bispectrum
>>> from astropy.io import fits
>>> moment0 = fits.open("Design4_21_0_0_flatrho_0021_13co.moment0.fits") # doctest: +SKIP
>>> bispec = Bispectrum(moment0) # doctest: +SKIP
>>> bispec.run(verbose=True, nsamples=1000) # doctest: +SKIP

Attributes Summary

bicoherence Bicoherence array.
bispectrum Bispectrum array.
bispectrum_logamp log amplitudes of the bispectrum.
data
distance
header
need_header_flag
no_data_flag
tracker Array showing the number of samples in each k_1 k_2 bin.

Methods Summary

azimuthal_slice(self, radii, delta_radii[, …]) Create an azimuthal slice of the bispectrum or bicoherence surfaces.
compute_bispectrum(self[, show_progress, …]) Do the computation.
input_data_header(self, data, header[, …]) Check if the header is given separately from the data type.
load_beam(self[, beam]) Try loading the beam from the header or a given object.
load_results(pickle_file) Load in a saved pickle file.
plot_surface(self[, save_name, show_bicoh, …]) Plot the bispectrum amplitude and (optionally) the bicoherence.
radial_slice(self, thetas, delta_thetas[, …]) Create a radial slice of the bispectrum (or bicoherence) plane.
run(self[, show_progress, use_pyfftw, …]) Compute the bispectrum.
save_results(self, output_name[, keep_data]) Save the results of the SCF to avoid re-computing.

Attributes Documentation

bicoherence

Bicoherence array.

bispectrum

Bispectrum array.

bispectrum_logamp

log amplitudes of the bispectrum.

data
distance
header
need_header_flag = True
no_data_flag = False
tracker

Array showing the number of samples in each k_1 k_2 bin.

Methods Documentation

azimuthal_slice(self, radii, delta_radii, bin_width=<Quantity 5. deg>, value='bispectrum', return_masks=False)[source] [edit on github]

Create an azimuthal slice of the bispectrum or bicoherence surfaces.

Parameters:
radii : float or np.ndarray

Radii in the bispectrum plane to extract slices at. Multiple slices are returned if radii is an array.

delta_radii : float or np.ndarray

The width around the radii in the bispectrum plane. If multiple radii are given, delta_radii must match the length of radii.

bin_width : Quantity, optional

The angular size of the bins used to create the slice.

value : str, optional

Which surface to create a profile from. Can be “bispectrum” (default), “bispectrum_logamp”, or “bicoherence”.

return_masks : bool, optional

Return the radial masks used to create the slices.

Returns:
azimuthal_slices : dict

Dictionary with the azimuthal slices. Each radius given in radii is a key in the dictionary. Each slice is a numpy array containing the bin centers, the averaged values, and the standard deviations in the azimuthal bins.

compute_bispectrum(self, show_progress=True, use_pyfftw=False, threads=1, nsamples=100, seed=1000, mean_subtract=False, **pyfftw_kwargs)[source] [edit on github]

Do the computation.

Parameters:
show_progress : optional, bool

Show progress bar while sampling the bispectrum.

use_pyfftw : bool, optional

Enable to use pyfftw, if it is installed.

threads : int, optional

Number of threads to use in FFT when using pyfftw.

nsamples : int, optional

Sets the number of samples to take at each vector magnitude.

seed : int, optional

Sets the seed for the distribution draws.

mean_subtract : bool, optional

Subtract the mean from the data before computing. This removes the “zero frequency” (i.e., constant) portion of the power, resulting in a loss of phase coherence along the k_1=k_2 line.

pyfft_kwargs : Passed to

rfft_to_fft. See here for a list of accepted kwargs.

input_data_header(self, data, header, need_copy=False) [edit on github]

Check if the header is given separately from the data type.

load_beam(self, beam=None) [edit on github]

Try loading the beam from the header or a given object.

Parameters:
beam : Beam, optional

The beam.

static load_results(pickle_file) [edit on github]

Load in a saved pickle file.

Parameters:
pickle_file : str

Name of filename to load in.

Returns:
self : Save statistic class

Statistic instance with saved results.

Examples

Load saved results. >>> stat = Statistic.load_results(“stat_saved.pkl”) # doctest: +SKIP

plot_surface(self, save_name=None, show_bicoh=True, cmap='viridis', contour_color='k')[source] [edit on github]

Plot the bispectrum amplitude and (optionally) the bicoherence.

Parameters:
save_name : str, optional

Save name for the figure. Enables saving the plot.

show_bicoh : bool, optional

Plot the bicoherence surface. Enabled by default.

cmap : {str, matplotlib color map}, optional

Colormap to use in the plots. Default is viridis.

contour_color : {str, RGB tuple}, optional

Color of the amplitude contours.

radial_slice(self, thetas, delta_thetas, bin_width=1.0, value='bispectrum', return_masks=False)[source] [edit on github]

Create a radial slice of the bispectrum (or bicoherence) plane.

Parameters:
thetas : Quantity

Azimuthal angles in the bispectrum plane to extract slices at. Multiple slices are returned if thetas is an array.

delta_thetas : Quantity

The width around the angle in the bispectrum plane. If multiple angles are given, delta_thetas must match the length of thetas.

bin_width : float, optional

The radial size of the bins used to create the slice.

value : str, optional

Which surface to create a profile from. Can be “bispectrum” (default), “bispectrum_logamp”, or “bicoherence”.

return_masks : bool, optional

Return the radial masks used to create the slices.

Returns:
radial_slices : dict

Dictionary with the radial slices. Each angle given in thetas is a key in the dictionary. Each slice is a numpy array containing the bin centers, the averaged values, and the standard deviations in the radial bins.

run(self, show_progress=True, use_pyfftw=False, threads=1, nsamples=100, seed=1000, mean_subtract=False, verbose=False, save_name=None, **pyfftw_kwargs)[source] [edit on github]

Compute the bispectrum. Necessary to maintain package standards.

Parameters:
show_progress : optional, bool

Show progress bar while sampling the bispectrum.

use_pyfftw : bool, optional

Enable to use pyfftw, if it is installed.

threads : int, optional

Number of threads to use in FFT when using pyfftw.

nsamples : int, optional

See compute_bispectrum.

seed : int, optional

See compute_bispectrum.

mean_subtract : bool, optional

See compute_bispectrum.

verbose : bool, optional

Enables plotting.

save_name : str,optional

Save the figure when a file name is given.

pyfftw_kwargs : Passed to

rfft_to_fft. See here for a list of accepted kwargs.

save_results(self, output_name, keep_data=False) [edit on github]

Save the results of the SCF to avoid re-computing. The pickled file will not include the data cube by default.

Parameters:
output_name : str

Name of the outputted pickle file.

keep_data : bool, optional

Save the data cube in the pickle file when enabled.