SCF

class turbustat.statistics.SCF(cube, header=None, size=11, roll_lags=None, distance=None)[source]

Bases: BaseStatisticMixIn

Computes the Spectral Correlation Function of a data cube (Rosolowsky et al, 1999).

Parameters:
cubenumpy.ndarray or astropy.io.fits.PrimaryHDU or astropy.io.fits.ImageHDU or SpectralCube

Data cube.

headerFITS header, optional

Header for the cube.

sizeint, optional

The total size of the lags used in one dimension in pixels. The maximum lag size will be (size - 1) / 2 in each direction.

roll_lagsndarray or Quantity, optional

Pass a custom array of lag values. An odd number of lags, centered at 0, must be given. If no units are given, it is assumed that the lags are in pixels. The lags should have symmetric positive and negative values (e.g., [-1, 0, 1]).

distanceQuantity, optional

Physical distance to the region in the data.

Examples

>>> from spectral_cube import SpectralCube
>>> from turbustat.statistics import SCF
>>> cube = SpectralCube.read("Design4.13co.fits")  
>>> scf = SCF(cube)  
>>> scf.run(verbose=True)  

Attributes Summary

data

distance

ellip2D

Fitted ellipticity of the 2D power-law.

ellip2D_err

Ellipticity standard error of the 2D power-law.

header

lags

Values of the lags, in pixels, to compute SCF at

need_header_flag

no_data_flag

roll_lags

Pixel values that the cube is rolled by to compute the SCF correlation surface.

scf_spectrum

Azimuthally averaged 1D SCF spectrum

scf_spectrum_stddev

Standard deviation of the scf_spectrum

scf_surface

SCF correlation array

slope

SCF spectrum slope

slope2D

Fitted slope of the 2D power-law.

slope2D_err

Slope standard error of the 2D power-law.

slope_err

1-sigma error on the SCF spectrum slope

theta2D

Fitted position angle of the 2D power-law.

theta2D_err

Position angle standard error of the 2D power-law.

xhigh

Upper limit for lags to consider in fits.

xlow

Lower limit for lags to consider in fits.

Methods Summary

compute_spectrum(**kwargs)

Compute the 1D spectrum as a function of lag.

compute_surface([boundary, show_progress])

Computes the SCF up to the given lag value.

fit_2Dplaw([fit_method, p0, xlow, xhigh, ...])

Model the 2D power-spectrum surface with an elliptical power-law model.

fit_plaw([xlow, xhigh, verbose, bootstrap])

Fit a power-law to the SCF spectrum.

fitted_model(xvals)

Computes the modelled power-law using the given x values.

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.

plot_fit([save_name, show_radial, ...])

Plot the SCF surface, radial profiles, and associated fits.

run([boundary, show_progress, xlow, xhigh, ...])

Computes all SCF outputs.

save_results(output_name[, keep_data])

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

Attributes Documentation

data
distance
ellip2D

Fitted ellipticity of the 2D power-law.

ellip2D_err

Ellipticity standard error of the 2D power-law.

header
lags

Values of the lags, in pixels, to compute SCF at

need_header_flag = True
no_data_flag = False
roll_lags

Pixel values that the cube is rolled by to compute the SCF correlation surface.

scf_spectrum

Azimuthally averaged 1D SCF spectrum

scf_spectrum_stddev

Standard deviation of the scf_spectrum

scf_surface

SCF correlation array

slope

SCF spectrum slope

slope2D

Fitted slope of the 2D power-law.

slope2D_err

Slope standard error of the 2D power-law.

slope_err

1-sigma error on the SCF spectrum slope

theta2D

Fitted position angle of the 2D power-law.

theta2D_err

Position angle standard error of the 2D power-law.

xhigh

Upper limit for lags to consider in fits.

xlow

Lower limit for lags to consider in fits.

Methods Documentation

compute_spectrum(**kwargs)[source]

Compute the 1D spectrum as a function of lag. Can optionally use log-spaced bins. kwargs are passed into the pspec function, which provides many options. The default settings are applicable in nearly all use cases.

Parameters:
kwargspassed to turbustat.statistics.psds.pspec
compute_surface(boundary='continuous', show_progress=True)[source]

Computes the SCF up to the given lag value. This is an expensive operation and could take a long time to calculate.

Parameters:
boundary{“continuous”, “cut”}

Treat the boundary as continuous (wrap-around) or cut values beyond the edge (i.e., for most observational data).

show_progressbool, optional

Show a progress bar when computing the surface. =

fit_2Dplaw(fit_method='LevMarq', p0=(), xlow=None, xhigh=None, bootstrap=True, niters=100, use_azimmask=False)[source]

Model the 2D power-spectrum surface with an elliptical power-law model.

Parameters:
fit_methodstr, optional

The algorithm fitting to use. Only ‘LevMarq’ is currently available.

p0tuple, optional

Initial parameters for fitting. If no values are given, the initial parameters start from the 1D fit parameters.

xlowQuantity, optional

Lower lag value limit to consider in the fit.

xhighQuantity, optional

Upper lag value limit to consider in the fit.

bootstrapbool, optional

Bootstrap using the model residuals to estimate the parameter standard errors. This tends to give more realistic intervals than the covariance matrix.

nitersint, optional

Number of bootstrap iterations.

use_azimmaskbool, optional

Use the azimuthal mask defined for the 1D spectrum, when azimuthal limit have been given.

fit_plaw(xlow=None, xhigh=None, verbose=False, bootstrap=False, **bootstrap_kwargs)[source]

Fit a power-law to the SCF spectrum.

Parameters:
xlowQuantity, optional

Lower lag value limit to consider in the fit.

xhighQuantity, optional

Upper lag value limit to consider in the fit.

verbosebool, optional

Show fit summary when enabled.

fitted_model(xvals)[source]

Computes the modelled power-law using the given x values.

Parameters:
xvalsQuantity

Values of lags to compute the model at.

Returns:
model_valuesndarray

Values of the model at the given values. Equivalent to log values of the SCF spectrum.

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

plot_fit(save_name=None, show_radial=True, show_residual=True, show_surface=True, contour_color='k', cmap='viridis', data_color='r', fit_color='k', xunit=Unit('pix'))[source]

Plot the SCF surface, radial profiles, and associated fits.

Parameters:
save_namestr, optional

Save name for the figure. Enables saving the plot.

show_radialbool, optional

Show the azimuthally-averaged 1D SCF spectrum and fit.

show_surfacebool, optional

Show the SCF surface and (if performed) fit.

show_residualbool, optional

Plot the residuals for the 1D SCF fit.

contour_color{str, RGB tuple}, optional

Color of the 2D fit contours.

cmap{str, matplotlib color map}, optional

Colormap to use in the plots. Default is viridis.

data_color{str, RGB tuple}, optional

Color of the azimuthally-averaged data.

fit_color{str, RGB tuple}, optional

Color of the 1D fit.

xunitUnit, optional

Choose the angular unit to convert to when ang_units is enabled.

run(boundary='continuous', show_progress=True, xlow=None, xhigh=None, fit_kwargs={}, fit_2D=True, fit_2D_kwargs={}, radialavg_kwargs={}, verbose=False, xunit=Unit('pix'), save_name=None)[source]

Computes all SCF outputs.

Parameters:
boundary{“continuous”, “cut”}

Treat the boundary as continuous (wrap-around) or cut values beyond the edge (i.e., for most observational data).

show_progressbool, optional

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

xlowQuantity, optional

See fit_plaw.

xhighQuantity, optional

See fit_plaw.

fit_kwargsdict, optional

Keyword arguments for SCF.fit_plaw. Use the xlow and xhigh keywords to provide fit limits.

fit_2Dbool, optional

Fit an elliptical power-law model to the 2D spectrum.

fit_2D_kwargsdict, optional

Keyword arguments for SCF.fit_2Dplaw. Use the xlow and xhigh keywords to provide fit limits.

radialavg_kwargsdict, optional

Passed to compute_spectrum.

verbosebool, optional

Enables plotting.

xunitUnit, optional

Choose the angular unit to convert to when ang_units is enabled.

save_namestr, optional

Save the figure when a file name is given.

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.