turbustat.statistics.
DeltaVariance
(img, header=None, weights=None, diam_ratio=1.5, lags=None, nlags=25)[source] [edit on github]¶Bases: turbustat.statistics.base_statistic.BaseStatisticMixIn
The delta-variance technique as described in Ossenkopf et al. (2008).
Parameters: | img : numpy.ndarray or astropy.io.fits.PrimaryHDU or spectral_cube.LowerDimensionalObject
header : FITS header, optional
weights : numpy.ndarray or astropy.io.fits.PrimaryHDU or spectral_cube.LowerDimensionalObject
diam_ratio : float, optional
lags : numpy.ndarray or list, optional
nlags : int, optional
|
---|
Attributes Summary
delta_var |
Delta Variance values. |
delta_var_error |
1-sigma errors on the Delta variance values. |
fit_range |
|
slope |
|
slope_err |
|
weights |
Array of weights. |
Methods Summary
compute_deltavar () |
Computes the delta-variance values and errors. |
do_convolutions ([allow_huge, boundary]) |
Perform the convolutions at all lags. |
fit_plaw ([xlow, xhigh, verbose]) |
Fit a power-law to the SCF spectrum. |
fitted_model (xvals) |
Computes the fitted power-law in log-log space using the given x values. |
run ([verbose, ang_units, unit, allow_huge, ...]) |
Compute the delta-variance. |
Attributes Documentation
delta_var
¶Delta Variance values.
delta_var_error
¶1-sigma errors on the Delta variance values.
fit_range
¶slope
¶slope_err
¶weights
¶Array of weights.
Methods Documentation
compute_deltavar
()[source] [edit on github]¶Computes the delta-variance values and errors.
do_convolutions
(allow_huge=False, boundary='wrap')[source] [edit on github]¶Perform the convolutions at all lags.
Parameters: | allow_huge : bool, optional
boundary : {“wrap”, “fill”}, optional
|
---|
fit_plaw
(xlow=None, xhigh=None, verbose=False)[source] [edit on github]¶Fit a power-law to the SCF spectrum.
Parameters: | xlow : float, optional
xhigh : float, optional
verbose : bool, optional
|
---|
fitted_model
(xvals)[source] [edit on github]¶Computes the fitted power-law in log-log space using the given x values.
Parameters: | xvals :
|
---|---|
Returns: | model_values :
|
run
(verbose=False, ang_units=False, unit=Unit("deg"), allow_huge=False, boundary='wrap', xlow=None, xhigh=None, save_name=None)[source] [edit on github]¶Compute the delta-variance.
Parameters: | verbose : bool, optional
ang_units : bool, optional
unit : u.Unit, optional
allow_huge : bool, optional
boundary : {“wrap”, “fill”}, optional
xlow : float, optional
xhigh : float, optional
save_name : str,optional
|
---|