turbustat.statistics.
PCA
(cube, n_eigs=None, distance=None)[source] [edit on github]¶Bases: turbustat.statistics.base_statistic.BaseStatisticMixIn
Implementation of Principal Component Analysis (Heyer & Brunt, 2002)
Parameters: | cube : numpy.ndarray or astropy.io.fits.PrimaryHDU or SpectralCube
n_eigs : int
distance :
|
---|
Examples
>>> from turbustat.statistics import PCA
>>> from spectral_cube import SpectralCube
>>> import astropy.units as u
>>> cube = SpectralCube.read("adv.fits")
>>> pca = PCA(cube, distance=250 * u.pc)
>>> pca.run(verbose=True)
Attributes Summary
eigvals |
All eigenvalues. |
eigvecs |
All eigenvectors. |
gamma |
Slope of the size-linewidth relation with correction from Chris Brunt’s thesis. |
gamma_error_range |
One-sigma error bounds on gamma. |
index |
Power-law index. |
index_error_range |
One-sigma error bounds on the index. |
intercept |
Intercept from the fits, converted out of the log value. |
intercept_error_range |
One-sigma error bounds on the intercept. |
n_eigs |
|
spatial_width |
Spatial widths for the first n_eigs components. |
spatial_width_error |
The 1-sigma error bounds on the spatial widths for the first n_eigs components. |
spectral_width |
Spectral widths for the first n_eigs components. |
spectral_width_error |
The error bounds on the spectral widths for the first n_eigs components. |
total_variance |
Total variance of all eigenvalues. |
var_proportion |
Proportion of variance described by the first n_eigs eigenvalues. |
Methods Summary
autocorr_images ([n_eigs]) |
Create the autocorrelation of the eigenimages. |
autocorr_spec ([n_eigs]) |
Create the autocorrelation spectra of the eigenvectors. |
compute_pca ([mean_sub, n_eigs, min_eigval, ...]) |
Create the covariance matrix and its eigenvalues. |
eigimages ([n_eigs]) |
Create eigenimages up to the n_eigs. |
find_spatial_widths ([method, ...]) |
Derive the spatial widths using the autocorrelation of the eigenimages. |
find_spectral_widths ([method, physical_units]) |
Derive the spectral widths using the autocorrelation of the eigenvectors. |
fit_plaw ([fit_method, verbose]) |
Fit the size-linewidth relation. |
model (x) |
Model with the fit parameters from fit_plaw |
noise_ACF ([n_eigs]) |
Create the noise autocorrelation function based off of the eigenvalues beyond PCA.n_eigs . |
run ([verbose, save_name, mean_sub, ...]) |
Run the decomposition and fitting in one step. |
sonic_length ([T_k, mu, use_gamma]) |
Estimate of the sonic length based on a given temperature. |
Attributes Documentation
eigvals
¶All eigenvalues.
eigvecs
¶All eigenvectors.
gamma_error_range
¶One-sigma error bounds on gamma.
index
¶Power-law index.
index_error_range
¶One-sigma error bounds on the index.
intercept
¶Intercept from the fits, converted out of the log value.
intercept_error_range
¶One-sigma error bounds on the intercept.
n_eigs
¶spatial_width_error
¶The 1-sigma error bounds on the spatial widths for the first
n_eigs
components.
total_variance
¶Total variance of all eigenvalues.
Methods Documentation
autocorr_images
(n_eigs=None)[source] [edit on github]¶Create the autocorrelation of the eigenimages.
Parameters: | n_eigs : None or int
|
---|---|
Returns: | acors : np.ndarray
|
autocorr_spec
(n_eigs=None)[source] [edit on github]¶Create the autocorrelation spectra of the eigenvectors.
Parameters: | n_eigs : None or int
|
---|---|
Returns: | acors : np.ndarray
|
compute_pca
(mean_sub=False, n_eigs='auto', min_eigval=None, eigen_cut_method='value')[source] [edit on github]¶Create the covariance matrix and its eigenvalues.
If mean_sub
is disabled, the first eigenvalue is dominated by the
mean of the data, not the variance.
Parameters: | mean_sub : bool, optional
n_eigs : {int, ‘auto’}, optional
min_eigval : float, optional
eigen_cut_method : {‘proportion’, ‘value’}, optional
|
---|
eigimages
(n_eigs=None)[source] [edit on github]¶Create eigenimages up to the n_eigs.
Parameters: | n_eigs : None or int
|
---|---|
Returns: | eigimgs :
|
find_spatial_widths
(method='contour', brunt_beamcorrect=True, beam_fwhm=None, physical_scales=True, distance=None)[source] [edit on github]¶Derive the spatial widths using the autocorrelation of the eigenimages.
Parameters: | methods : {‘contour’, ‘fit’, ‘interpolate’, ‘xinterpolate’}, optional
brunt_beamcorrect : bool, optional
beam_fwhm : None of
physical_scales : bool, optional
distance :
|
---|
find_spectral_widths
(method='walk-down', physical_units=True)[source] [edit on github]¶Derive the spectral widths using the autocorrelation of the eigenvectors.
Parameters: | method : {“walk-down”, “fit”, “interpolate”}, optional
physical_units : bool, optional
|
---|
fit_plaw
(fit_method='odr', verbose=False, **kwargs)[source] [edit on github]¶Fit the size-linewidth relation. This is done through Orthogonal Distance Regression (via scipy), or through MCMC (requires installing emcee).
Parameters: | fit_method : {‘odr’, ‘bayes’}, optional
verbose : bool, optional
**kwargs
|
---|
model
(x)[source] [edit on github]¶Model with the fit parameters from fit_plaw
noise_ACF
(n_eigs=-10)[source] [edit on github]¶Create the noise autocorrelation function based off of the eigenvalues
beyond PCA.n_eigs
.
run
(verbose=False, save_name=None, mean_sub=False, decomp_only=False, n_eigs='auto', min_eigval=None, eigen_cut_method='value', spatial_method='contour', spectral_method='walk-down', fit_method='odr', beam_fwhm=None, brunt_beamcorrect=True)[source] [edit on github]¶Run the decomposition and fitting in one step.
Parameters: | verbose : bool, optional
save_name : str,optional
mean_sub : bool, optional
decomp_only : bool, optional
n_eigs : {“auto”, int}, optional
min_eigval : float, optional
eigen_cut_method : {‘proportion’, ‘value’}, optional
spatial_method : str, optional
spectral_method : str, optional
fit_method : str, optional
beam_fwhm : None of
brunt_beamcorrect : bool, optional
|
---|
sonic_length
(T_k=<Quantity 10.0 K>, mu=1.36, use_gamma=True)[source] [edit on github]¶Estimate of the sonic length based on a given temperature. Uses the intercept from the fit.
Based on sonic.pro used in the Heyer & Brunt PCA implementation.
Because error from the MCMC fit need not be symmetric, the MCMC samples are needed to provide the correct CIs for the sonic length.
Parameters: | T_k :
mu : float, optional
use_gamma : bool, optional |
---|---|
Returns: | lambd :
lambd_error_range :
|