StatMoments

class turbustat.statistics.StatMoments(img, header=None, weights=None, radius=<Quantity 5. pix>, nbins=None, distance=None)[source]

Bases: BaseStatisticMixIn

Statistical Moments of an image. See Burkhart et al. (2010) for the methods used. By specifying the radius of circular mask, the mean, variance, skewness, and kurtosis are calculated within the circular mask for every pixel in the image. The distributions of these moments can be compared between data sets.

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

2D image.

headerFITS header, optional

The image header. Needed for the pixel scale.

weightsnumpy.ndarray or astropy.io.fits.PrimaryHDU or astropy.io.fits.ImageHDU or spectral_cube.Projection or spectral_cube.Slice

2D array of weights. Uniform weights are used if none are given.

radiusQuantity, optional

Radius of circle to use when computing moments. When angular or physical units are given, they will be rounded down to the nearest pixel size.

nbinsarray or int, optional

Number of bins to use in the histogram.

distanceQuantity, optional

Physical distance to the region in the data.

Attributes Summary

data

distance

header

kurtosis_array

The array of local kurtosiss.

kurtosis_extrema

The extrema of the kurtosis array.

kurtosis_hist

The histogram bins and values for the kurtosis array.

mean_array

The array of local means.

mean_extrema

The extrema of the mean array.

mean_hist

The histogram bins and values for the mean array.

need_header_flag

no_data_flag

radius

skewness_array

The array of local skewnesss.

skewness_extrema

The extrema of the skewness array.

skewness_hist

The histogram bins and values for the skewness array.

variance_array

The array of local variances.

variance_extrema

The extrema of the variance array.

variance_hist

The histogram bins and values for the variance array.

Methods Summary

array_moments()

Moments over the entire image.

compute_spatial_distrib([radius, periodic, ...])

Compute the moments over circular region with the specified radius.

input_data_header(data, header[, need_copy])

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

load_beam([beam])

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

load_results(pickle_file)

Load in a saved pickle file.

make_spatial_histograms([mean_bins, ...])

Create histograms of the moments.

plot_histograms([new_figure, save_name, ...])

Plot the histograms of each moment.

plot_maps([save_name, cmap, contour_cmap])

Plot the maps of locally-estimated moments.

run([show_progress, verbose, save_name, ...])

Compute the entire method.

save_results(output_name[, keep_data])

Save the results of the SCF to avoid re-computing.

Attributes Documentation

data
distance
header
kurtosis_array

The array of local kurtosiss.

kurtosis_extrema

The extrema of the kurtosis array.

kurtosis_hist

The histogram bins and values for the kurtosis array. The first element is the array of bins, and the second contains the values.

mean_array

The array of local means.

mean_extrema

The extrema of the mean array.

mean_hist

The histogram bins and values for the mean array. The first element is the array of bins, and the second contains the values.

need_header_flag = True
no_data_flag = False
radius
skewness_array

The array of local skewnesss.

skewness_extrema

The extrema of the skewness array.

skewness_hist

The histogram bins and values for the skewness array. The first element is the array of bins, and the second contains the values.

variance_array

The array of local variances.

variance_extrema

The extrema of the variance array.

variance_hist

The histogram bins and values for the variance array. The first element is the array of bins, and the second contains the values.

Methods Documentation

array_moments()[source]

Moments over the entire image.

compute_spatial_distrib(radius=None, periodic=True, min_frac=0.8, show_progress=True)[source]

Compute the moments over circular region with the specified radius.

Parameters:
radiusQuantity, optional

Override the radius size of the region.

periodicbool, optional

Specify whether the boundaries can be wrapped. Default is True.

min_fracfloat, optional

A number between 0 and 1 that sets the minimum fraction of data in each region that are finite. A value of 1.0 requires that no NaNs be in the region.

show_progressbool, optional

Show a progress bar during the creation of the covariance matrix.

input_data_header(data, header, need_copy=False)

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

load_beam(beam=None)

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

Parameters:
beamBeam, optional

The beam.

static load_results(pickle_file)

Load in a saved pickle file.

Parameters:
pickle_filestr

Name of filename to load in.

Returns:
selfSave statistic class

Statistic instance with saved results.

Examples

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

make_spatial_histograms(mean_bins=None, variance_bins=None, skewness_bins=None, kurtosis_bins=None)[source]

Create histograms of the moments. If an optional set of bins is not given, \(\sqrt{N}\) equally-size bins will be created, where \(N\) is the number of elements in the array. The histogram values are normalized so that the sum of the values in the bins, multiplied by the bin width is 1.

Parameters:
mean_binsarray, optional

Bins to use for the histogram of the mean array.

variance_binsarray, optional

Bins to use for the histogram of the variance array.

skewness_binsarray, optional

Bins to use for the histogram of the skewness array.

kurtosis_binsarray, optional

Bins to use for the histogram of the kurtosis array.

plot_histograms(new_figure=True, save_name=None, hist_color='r', face_color='k')[source]

Plot the histograms of each moment.

Parameters:
new_figurebool, optional

Creates a new matplotlib figure.

save_namestr, optional

The filename to save the plot as. This enables saving of the plot.

plot_maps(save_name=None, cmap='binary', contour_cmap='viridis')[source]

Plot the maps of locally-estimated moments.

Parameters:
save_namestr, optional

Save name for the figure. Enables saving the plot.

cmap{str, matplotlib colormap}, optional

Colormap for the images.

contour_cmap{str, matplotlib colormap}, optional

Colormap for the contours.

run(show_progress=True, verbose=False, save_name=None, radius=None, periodic=True, min_frac=0.8, **hist_kwargs)[source]

Compute the entire method.

Parameters:
show_progressbool, optional

Show a progress bar during the creation of the covariance matrix.

verbosebool, optional

Enables plotting.

save_namestr,optional

Save the figure when a file name is given.

radiusQuantity, optional

Override the radius size of the region.

periodicbool, optional

Specify whether the boundaries can be wrapped. Default is True.

min_fracfloat, optional

A number between 0 and 1 that sets the minimum fraction of data in each region that are finite. A value of 1.0 requires that no NaNs be in the region.

hist_kwargsPassed to make_spatial_histograms.
save_results(output_name, keep_data=False)

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

Parameters:
output_namestr

Name of the outputted pickle file.

keep_databool, optional

Save the data cube in the pickle file when enabled.