LombScargle¶
-
class
astropy.timeseries.
LombScargle
(t, y, dy=None, fit_mean=True, center_data=True, nterms=1, normalization='standard')[source]¶ Bases:
astropy.timeseries.BasePeriodogram
Compute the Lomb-Scargle Periodogram.
This implementations here are based on code presented in [1] and [2]; if you use this functionality in an academic application, citation of those works would be appreciated.
- Parameters
- tarray_like or Quantity
sequence of observation times
- yarray_like or Quantity
sequence of observations associated with times t
- dyfloat, array_like, or Quantity, optional
error or sequence of observational errors associated with times t
- fit_meanbool, optional
if True, include a constant offset as part of the model at each frequency. This can lead to more accurate results, especially in the case of incomplete phase coverage.
- center_databool, optional
if True, pre-center the data by subtracting the weighted mean of the input data. This is especially important if fit_mean = False
- ntermsint, optional
number of terms to use in the Fourier fit
- normalization{‘standard’, ‘model’, ‘log’, ‘psd’}, optional
Normalization to use for the periodogram.
References
- 1
Vanderplas, J., Connolly, A. Ivezic, Z. & Gray, A. Introduction to astroML: Machine learning for astrophysics. Proceedings of the Conference on Intelligent Data Understanding (2012)
- 2
VanderPlas, J. & Ivezic, Z. Periodograms for Multiband Astronomical Time Series. ApJ 812.1:18 (2015)
Examples
Generate noisy periodic data:
>>> rand = np.random.RandomState(42) >>> t = 100 * rand.rand(100) >>> y = np.sin(2 * np.pi * t) + rand.randn(100)
Compute the Lomb-Scargle periodogram on an automatically-determined frequency grid & find the frequency of max power:
>>> frequency, power = LombScargle(t, y).autopower() >>> frequency[np.argmax(power)] 1.0016662310392956
Compute the Lomb-Scargle periodogram at a user-specified frequency grid:
>>> freq = np.arange(0.8, 1.3, 0.1) >>> LombScargle(t, y).power(freq) array([0.0204304 , 0.01393845, 0.35552682, 0.01358029, 0.03083737])
If the inputs are astropy Quantities with units, the units will be validated and the outputs will also be Quantities with appropriate units:
>>> from astropy import units as u >>> t = t * u.s >>> y = y * u.mag >>> frequency, power = LombScargle(t, y).autopower() >>> frequency.unit Unit("1 / s") >>> power.unit Unit(dimensionless)
Note here that the Lomb-Scargle power is always a unitless quantity, because it is related to the \(\chi^2\) of the best-fit periodic model at each frequency.
Attributes Summary
Methods Summary
autofrequency
(self[, samples_per_peak, …])Determine a suitable frequency grid for data.
autopower
(self[, method, method_kwds, …])Compute Lomb-Scargle power at automatically-determined frequencies.
design_matrix
(self, frequency[, t])Compute the design matrix for a given frequency
distribution
(self, power[, cumulative])Expected periodogram distribution under the null hypothesis.
false_alarm_level
(self, false_alarm_probability)Level of maximum at a given false alarm probability.
false_alarm_probability
(self, power[, …])False alarm probability of periodogram maxima under the null hypothesis.
from_timeseries
(timeseries[, …])Initialize a periodogram from a time series object.
model
(self, t, frequency)Compute the Lomb-Scargle model at the given frequency.
model_parameters
(self, frequency[, units])Compute the best-fit model parameters at the given frequency.
offset
(self)Return the offset of the model
power
(self, frequency[, normalization, …])Compute the Lomb-Scargle power at the given frequencies.
Attributes Documentation
-
available_methods
= ['auto', 'slow', 'chi2', 'cython', 'fast', 'fastchi2', 'scipy']¶
Methods Documentation
-
autofrequency
(self, samples_per_peak=5, nyquist_factor=5, minimum_frequency=None, maximum_frequency=None, return_freq_limits=False)[source]¶ Determine a suitable frequency grid for data.
Note that this assumes the peak width is driven by the observational baseline, which is generally a good assumption when the baseline is much larger than the oscillation period. If you are searching for periods longer than the baseline of your observations, this may not perform well.
Even with a large baseline, be aware that the maximum frequency returned is based on the concept of “average Nyquist frequency”, which may not be useful for irregularly-sampled data. The maximum frequency can be adjusted via the nyquist_factor argument, or through the maximum_frequency argument.
- Parameters
- samples_per_peakfloat, optional
The approximate number of desired samples across the typical peak
- nyquist_factorfloat, optional
The multiple of the average nyquist frequency used to choose the maximum frequency if maximum_frequency is not provided.
- minimum_frequencyfloat, optional
If specified, then use this minimum frequency rather than one chosen based on the size of the baseline.
- maximum_frequencyfloat, optional
If specified, then use this maximum frequency rather than one chosen based on the average nyquist frequency.
- return_freq_limitsbool, optional
if True, return only the frequency limits rather than the full frequency grid.
- Returns
- frequencyndarray or Quantity
The heuristically-determined optimal frequency bin
-
autopower
(self, method='auto', method_kwds=None, normalization=None, samples_per_peak=5, nyquist_factor=5, minimum_frequency=None, maximum_frequency=None)[source]¶ Compute Lomb-Scargle power at automatically-determined frequencies.
- Parameters
- methodstr, optional
specify the lomb scargle implementation to use. Options are:
‘auto’: choose the best method based on the input
‘fast’: use the O[N log N] fast method. Note that this requires evenly-spaced frequencies: by default this will be checked unless
assume_regular_frequency
is set to True.‘slow’: use the O[N^2] pure-python implementation
‘cython’: use the O[N^2] cython implementation. This is slightly faster than method=’slow’, but much more memory efficient.
‘chi2’: use the O[N^2] chi2/linear-fitting implementation
‘fastchi2’: use the O[N log N] chi2 implementation. Note that this requires evenly-spaced frequencies: by default this will be checked unless
assume_regular_frequency
is set to True.‘scipy’: use
scipy.signal.lombscargle
, which is an O[N^2] implementation written in C. Note that this does not support heteroskedastic errors.
- method_kwdsdict, optional
additional keywords to pass to the lomb-scargle method
- normalization{‘standard’, ‘model’, ‘log’, ‘psd’}, optional
If specified, override the normalization specified at instantiation.
- samples_per_peakfloat, optional
The approximate number of desired samples across the typical peak
- nyquist_factorfloat, optional
The multiple of the average nyquist frequency used to choose the maximum frequency if maximum_frequency is not provided.
- minimum_frequencyfloat, optional
If specified, then use this minimum frequency rather than one chosen based on the size of the baseline.
- maximum_frequencyfloat, optional
If specified, then use this maximum frequency rather than one chosen based on the average nyquist frequency.
- Returns
- frequency, powerndarrays
The frequency and Lomb-Scargle power
-
design_matrix
(self, frequency, t=None)[source]¶ Compute the design matrix for a given frequency
- Parameters
- Returns
- Xnp.ndarray (len(t), n_parameters)
The design matrix for the model at the given frequency.
See also
-
distribution
(self, power, cumulative=False)[source]¶ Expected periodogram distribution under the null hypothesis.
This computes the expected probability distribution or cumulative probability distribution of periodogram power, under the null hypothesis of a non-varying signal with Gaussian noise. Note that this is not the same as the expected distribution of peak values; for that see the
false_alarm_probability()
method.- Parameters
- powerarray_like
The periodogram power at which to compute the distribution.
- cumulativebool, optional
If True, then return the cumulative distribution.
- Returns
- distnp.ndarray
The probability density or cumulative probability associated with the provided powers.
See also
-
false_alarm_level
(self, false_alarm_probability, method='baluev', samples_per_peak=5, nyquist_factor=5, minimum_frequency=None, maximum_frequency=None, method_kwds=None)[source]¶ Level of maximum at a given false alarm probability.
This gives an estimate of the periodogram level corresponding to a specified false alarm probability for the largest peak, assuming a null hypothesis of non-varying data with Gaussian noise.
- Parameters
- false_alarm_probabilityarray_like
The false alarm probability (0 < fap < 1).
- maximum_frequencyfloat
The maximum frequency of the periodogram.
- method{‘baluev’, ‘davies’, ‘naive’, ‘bootstrap’}, optional
The approximation method to use; default=’baluev’.
- method_kwdsdict, optional
Additional method-specific keywords.
- Returns
- powernp.ndarray
The periodogram peak height corresponding to the specified false alarm probability.
See also
Notes
The true probability distribution for the largest peak cannot be determined analytically, so each method here provides an approximation to the value. The available methods are:
“baluev” (default): the upper-limit to the alias-free probability, using the approach of Baluev (2008) [1].
“davies” : the Davies upper bound from Baluev (2008) [1].
“naive” : the approximate probability based on an estimated effective number of independent frequencies.
“bootstrap” : the approximate probability based on bootstrap resamplings of the input data.
Note also that for normalization=’psd’, the distribution can only be computed for periodograms constructed with errors specified.
References
-
false_alarm_probability
(self, power, method='baluev', samples_per_peak=5, nyquist_factor=5, minimum_frequency=None, maximum_frequency=None, method_kwds=None)[source]¶ False alarm probability of periodogram maxima under the null hypothesis.
This gives an estimate of the false alarm probability given the height of the largest peak in the periodogram, based on the null hypothesis of non-varying data with Gaussian noise.
- Parameters
- powerarray_like
The periodogram value.
- method{‘baluev’, ‘davies’, ‘naive’, ‘bootstrap’}, optional
The approximation method to use.
- maximum_frequencyfloat
The maximum frequency of the periodogram.
- method_kwdsdict, optional
Additional method-specific keywords.
- Returns
- false_alarm_probabilitynp.ndarray
The false alarm probability
See also
Notes
The true probability distribution for the largest peak cannot be determined analytically, so each method here provides an approximation to the value. The available methods are:
“baluev” (default): the upper-limit to the alias-free probability, using the approach of Baluev (2008) [1].
“davies” : the Davies upper bound from Baluev (2008) [1].
“naive” : the approximate probability based on an estimated effective number of independent frequencies.
“bootstrap” : the approximate probability based on bootstrap resamplings of the input data.
Note also that for normalization=’psd’, the distribution can only be computed for periodograms constructed with errors specified.
References
-
classmethod
from_timeseries
(timeseries, signal_column_name=None, uncertainty=None, **kwargs)¶ Initialize a periodogram from a time series object.
If a binned time series is passed, the time at the center of the bins is used. Also note that this method automatically gets rid of NaN/undefined values when initalizing the periodogram.
- Parameters
- signal_column_namestr
The name of the column containing the signal values to use.
- uncertaintystr or float or
Quantity
, optional The name of the column containing the errors on the signal, or the value to use for the error, if a scalar.
- **kwargs
Additional keyword arguments are passed to the initializer for this periodogram class.
-
model
(self, t, frequency)[source]¶ Compute the Lomb-Scargle model at the given frequency.
The model at a particular frequency is a linear model: model = offset + dot(design_matrix, model_parameters)
- Parameters
- tarray_like or Quantity, length n_samples
times at which to compute the model
- frequencyfloat
the frequency for the model
- Returns
- ynp.ndarray, length n_samples
The model fit corresponding to the input times
See also
-
model_parameters
(self, frequency, units=True)[source]¶ Compute the best-fit model parameters at the given frequency.
The model described by these parameters is:
\[y(t; f, \vec{\theta}) = \theta_0 + \sum_{n=1}^{\tt nterms} [\theta_{2n-1}\sin(2\pi n f t) + \theta_{2n}\cos(2\pi n f t)]\]where \(\vec{\theta}\) is the array of parameters returned by this function.
- Parameters
- frequencyfloat
the frequency for the model
- unitsbool
If True (default), return design matrix with data units.
- Returns
- thetanp.ndarray (n_parameters,)
The best-fit model parameters at the given frequency.
See also
-
offset
(self)[source]¶ Return the offset of the model
The offset of the model is the (weighted) mean of the y values. Note that if self.center_data is False, the offset is 0 by definition.
- Returns
- offsetscalar
See also
-
power
(self, frequency, normalization=None, method='auto', assume_regular_frequency=False, method_kwds=None)[source]¶ Compute the Lomb-Scargle power at the given frequencies.
- Parameters
- frequencyarray_like or Quantity
frequencies (not angular frequencies) at which to evaluate the periodogram. Note that in order to use method=’fast’, frequencies must be regularly-spaced.
- methodstr, optional
specify the lomb scargle implementation to use. Options are:
‘auto’: choose the best method based on the input
‘fast’: use the O[N log N] fast method. Note that this requires evenly-spaced frequencies: by default this will be checked unless
assume_regular_frequency
is set to True.‘slow’: use the O[N^2] pure-python implementation
‘cython’: use the O[N^2] cython implementation. This is slightly faster than method=’slow’, but much more memory efficient.
‘chi2’: use the O[N^2] chi2/linear-fitting implementation
‘fastchi2’: use the O[N log N] chi2 implementation. Note that this requires evenly-spaced frequencies: by default this will be checked unless
assume_regular_frequency
is set to True.‘scipy’: use
scipy.signal.lombscargle
, which is an O[N^2] implementation written in C. Note that this does not support heteroskedastic errors.
- assume_regular_frequencybool, optional
if True, assume that the input frequency is of the form freq = f0 + df * np.arange(N). Only referenced if method is ‘auto’ or ‘fast’.
- normalization{‘standard’, ‘model’, ‘log’, ‘psd’}, optional
If specified, override the normalization specified at instantiation.
- fit_meanbool, optional
If True, include a constant offset as part of the model at each frequency. This can lead to more accurate results, especially in the case of incomplete phase coverage.
- center_databool, optional
If True, pre-center the data by subtracting the weighted mean of the input data. This is especially important if fit_mean = False.
- method_kwdsdict, optional
additional keywords to pass to the lomb-scargle method
- Returns
- powerndarray
The Lomb-Scargle power at the specified frequency