SigmaClip¶
-
class
astropy.stats.SigmaClip(sigma=3.0, sigma_lower=None, sigma_upper=None, maxiters=5, cenfunc='median', stdfunc='std')[source]¶ Bases:
objectClass to perform sigma clipping.
The data will be iterated over, each time rejecting values that are less or more than a specified number of standard deviations from a center value.
Clipped (rejected) pixels are those where:
data < cenfunc(data [,axis=int]) - (sigma_lower * stdfunc(data [,axis=int])) data > cenfunc(data [,axis=int]) + (sigma_upper * stdfunc(data [,axis=int]))
Invalid data values (i.e. NaN or inf) are automatically clipped.
For a functional interface to sigma clipping, see
sigma_clip().Note
scipy.stats.sigmaclip provides a subset of the functionality in this class. Also, its input data cannot be a masked array and it does not handle data that contains invalid values (i.e. NaN or inf). Also note that it uses the mean as the centering function.
If your data is a
ndarraywith no invalid values and you want to use the mean as the centering function withaxis=Noneand iterate to convergence, thenscipy.stats.sigmaclipis ~25-30% faster than the equivalent settings here (s = SigmaClip(cenfunc='mean', maxiters=None); s(data, axis=None)).Parameters: - sigma : float, optional
The number of standard deviations to use for both the lower and upper clipping limit. These limits are overridden by
sigma_lowerandsigma_upper, if input. The default is 3.- sigma_lower : float or
None, optional The number of standard deviations to use as the lower bound for the clipping limit. If
Nonethen the value ofsigmais used. The default isNone.- sigma_upper : float or
None, optional The number of standard deviations to use as the upper bound for the clipping limit. If
Nonethen the value ofsigmais used. The default isNone.- maxiters : int or
None, optional The maximum number of sigma-clipping iterations to perform or
Noneto clip until convergence is achieved (i.e., iterate until the last iteration clips nothing). If convergence is achieved prior tomaxitersiterations, the clipping iterations will stop. The default is 5.- cenfunc : {‘median’, ‘mean’} or callable, optional
The statistic or callable function/object used to compute the center value for the clipping. If set to
'median'or'mean'then having the optional bottleneck package installed will result in the best performance. If using a callable function/object and theaxiskeyword is used, then it must be callable that can ignore NaNs (e.g.numpy.nanmean) and has anaxiskeyword to return an array with axis dimension(s) removed. The default is'median'.- stdfunc : {‘std’} or callable, optional
The statistic or callable function/object used to compute the standard deviation about the center value. If set to
'std'then having the optional bottleneck package installed will result in the best performance. If using a callable function/object and theaxiskeyword is used, then it must be callable that can ignore NaNs (e.g.numpy.nanstd) and has anaxiskeyword to return an array with axis dimension(s) removed. The default is'std'.
See also
Examples
This example uses a data array of random variates from a Gaussian distribution. We clip all points that are more than 2 sample standard deviations from the median. The result is a masked array, where the mask is
Truefor clipped data:>>> from astropy.stats import SigmaClip >>> from numpy.random import randn >>> randvar = randn(10000) >>> sigclip = SigmaClip(sigma=2, maxiters=5) >>> filtered_data = sigclip(randvar)
This example clips all points that are more than 3 sigma relative to the sample mean, clips until convergence, returns an unmasked
ndarray, and modifies the data in-place:>>> from astropy.stats import SigmaClip >>> from numpy.random import randn >>> from numpy import mean >>> randvar = randn(10000) >>> sigclip = SigmaClip(sigma=3, maxiters=None, cenfunc='mean') >>> filtered_data = sigclip(randvar, masked=False, copy=False)
This example sigma clips along one axis:
>>> from astropy.stats import SigmaClip >>> from numpy.random import normal >>> from numpy import arange, diag, ones >>> data = arange(5) + normal(0., 0.05, (5, 5)) + diag(ones(5)) >>> sigclip = SigmaClip(sigma=2.3) >>> filtered_data = sigclip(data, axis=0)
Note that along the other axis, no points would be clipped, as the standard deviation is higher.
Methods Summary
__call__(data[, axis, masked, …])Perform sigma clipping on the provided data. Methods Documentation
-
__call__(data, axis=None, masked=True, return_bounds=False, copy=True)[source]¶ Perform sigma clipping on the provided data.
Parameters: - data : array-like or
MaskedArray The data to be sigma clipped.
- axis :
Noneor int or tuple of int, optional The axis or axes along which to sigma clip the data. If
None, then the flattened data will be used.axisis passed to thecenfuncandstdfunc. The default isNone.- masked : bool, optional
If
True, then aMaskedArrayis returned, where the mask isTruefor clipped values. IfFalse, then andarrayand the minimum and maximum clipping thresholds are returned. The default isTrue.- return_bounds : bool, optional
If
True, then the minimum and maximum clipping bounds are also returned.- copy : bool, optional
If
True, then thedataarray will be copied. IfFalseandmasked=True, then the returned masked array data will contain the same array as the inputdata(ifdatais andarrayorMaskedArray). The default isTrue.
Returns: - result : flexible
If
masked=True, then aMaskedArrayis returned, where the mask isTruefor clipped values. Ifmasked=False, then andarrayis returned.If
return_bounds=True, then in addition to the (masked) array above, the minimum and maximum clipping bounds are returned.If
masked=Falseandaxis=None, then the output array is a flattened 1Dndarraywhere the clipped values have been removed. Ifreturn_bounds=Truethen the returned minimum and maximum thresholds are scalars.If
masked=Falseandaxisis specified, then the outputndarraywill have the same shape as the inputdataand containnp.nanwhere values were clipped. Ifreturn_bounds=Truethen the returned minimum and maximum clipping thresholds will be bendarrays.
- data : array-like or