LombScargle

class astropy.stats.LombScargle(t, y, dy=None, fit_mean=True, center_data=True, nterms=1, normalization='standard')[source] [edit on github]

Bases: object

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:
t : array_like or Quantity

sequence of observation times

y : array_like or Quantity

sequence of observations associated with times t

dy : float, array_like or Quantity (optional)

error or sequence of observational errors associated with times t

fit_mean : bool (optional, default=True)

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_data : bool (optional, default=True)

if True, pre-center the data by subtracting the weighted mean of the input data. This is especially important if fit_mean = False

nterms : int (optional, default=1)

number of terms to use in the Fourier fit

normalization : {‘standard’, ‘model’, ‘log’, ‘psd’}, optional

Normalization to use for the periodogram.

References

[1](1, 2) 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](1, 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)]  # doctest: +FLOAT_CMP
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)  # doctest: +FLOAT_CMP
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

available_methods

Methods Summary

autofrequency([samples_per_peak, …]) Determine a suitable frequency grid for data.
autopower([method, method_kwds, …]) Compute Lomb-Scargle power at automatically-determined frequencies.
distribution(power[, cumulative]) Expected periodogram distribution under the null hypothesis.
false_alarm_level(false_alarm_probability[, …]) Level of maximum at a given false alarm probability.
false_alarm_probability(power[, method, …]) False alarm probability of periodogram maxima under the null hypothesis.
model(t, frequency) Compute the Lomb-Scargle model at the given frequency.
power(frequency[, normalization, method, …]) Compute the Lomb-Scargle power at the given frequencies.

Attributes Documentation

available_methods = ['auto', 'slow', 'chi2', 'cython', 'fast', 'fastchi2', 'scipy']

Methods Documentation

autofrequency(samples_per_peak=5, nyquist_factor=5, minimum_frequency=None, maximum_frequency=None, return_freq_limits=False)[source] [edit on github]

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_peak : float (optional, default=5)

The approximate number of desired samples across the typical peak

nyquist_factor : float (optional, default=5)

The multiple of the average nyquist frequency used to choose the maximum frequency if maximum_frequency is not provided.

minimum_frequency : float (optional)

If specified, then use this minimum frequency rather than one chosen based on the size of the baseline.

maximum_frequency : float (optional)

If specified, then use this maximum frequency rather than one chosen based on the average nyquist frequency.

return_freq_limits : bool (optional)

if True, return only the frequency limits rather than the full frequency grid.

Returns:
frequency : ndarray or Quantity

The heuristically-determined optimal frequency bin

autopower(method='auto', method_kwds=None, normalization=None, samples_per_peak=5, nyquist_factor=5, minimum_frequency=None, maximum_frequency=None)[source] [edit on github]

Compute Lomb-Scargle power at automatically-determined frequencies.

Parameters:
method : string (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_kwds : dict (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_peak : float (optional, default=5)

The approximate number of desired samples across the typical peak

nyquist_factor : float (optional, default=5)

The multiple of the average nyquist frequency used to choose the maximum frequency if maximum_frequency is not provided.

minimum_frequency : float (optional)

If specified, then use this minimum frequency rather than one chosen based on the size of the baseline.

maximum_frequency : float (optional)

If specified, then use this maximum frequency rather than one chosen based on the average nyquist frequency.

Returns:
frequency, power : ndarrays

The frequency and Lomb-Scargle power

distribution(power, cumulative=False)[source] [edit on github]

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:
power : array_like

The periodogram power at which to compute the distribution.

cumulative : bool (optional)

If True, then return the cumulative distribution.

Returns:
dist : np.ndarray

The probability density or cumulative probability associated with the provided powers.

false_alarm_level(false_alarm_probability, method='baluev', samples_per_peak=5, nyquist_factor=5, minimum_frequency=None, maximum_frequency=None, method_kwds=None)[source] [edit on github]

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_probability : array-like

The false alarm probability (0 < fap < 1).

maximum_frequency : float

The maximum frequency of the periodogram.

method : {‘baluev’, ‘davies’, ‘naive’, ‘bootstrap’}, optional

The approximation method to use; default=’baluev’.

method_kwds : dict, optional

Additional method-specific keywords.

Returns:
power : np.ndarray

The periodogram peak height corresponding to the specified false alarm probability.

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

[1](1, 2, 3) Baluev, R.V. MNRAS 385, 1279 (2008)
false_alarm_probability(power, method='baluev', samples_per_peak=5, nyquist_factor=5, minimum_frequency=None, maximum_frequency=None, method_kwds=None)[source] [edit on github]

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:
power : array-like

The periodogram value.

method : {‘baluev’, ‘davies’, ‘naive’, ‘bootstrap’}, optional

The approximation method to use.

maximum_frequency : float

The maximum frequency of the periodogram.

method_kwds : dict (optional)

Additional method-specific keywords.

Returns:
false_alarm_probability : np.ndarray

The false alarm probability

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

[1](1, 2, 3) Baluev, R.V. MNRAS 385, 1279 (2008)
model(t, frequency)[source] [edit on github]

Compute the Lomb-Scargle model at the given frequency.

Parameters:
t : array_like or Quantity, length n_samples

times at which to compute the model

frequency : float

the frequency for the model

Returns:
y : np.ndarray, length n_samples

The model fit corresponding to the input times

power(frequency, normalization=None, method='auto', assume_regular_frequency=False, method_kwds=None)[source] [edit on github]

Compute the Lomb-Scargle power at the given frequencies.

Parameters:
frequency : array_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.

method : string (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_frequency : bool (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_mean : bool (optional, default=True)

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_data : bool (optional, default=True)

If True, pre-center the data by subtracting the weighted mean of the input data. This is especially important if fit_mean = False.

method_kwds : dict (optional)

additional keywords to pass to the lomb-scargle method

Returns:
power : ndarray

The Lomb-Scargle power at the specified frequency