BoxLeastSquares#

class astropy.timeseries.BoxLeastSquares(t, y, dy=None)[source]#

Bases: BasePeriodogram

Compute the box least squares periodogram.

This method is a commonly used tool for discovering transiting exoplanets or eclipsing binaries in photometric time series datasets. This implementation is based on the “box least squares (BLS)” method described in [1] and [2].

Parameters:
tarray_like, Quantity, Time, or TimeDelta

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.

References

[1]

Kovacs, Zucker, & Mazeh (2002), A&A, 391, 369 (arXiv:astro-ph/0206099)

[2]

Hartman & Bakos (2016), Astronomy & Computing, 17, 1 (arXiv:1605.06811)

Examples

Generate noisy data with a transit:

>>> rand = np.random.default_rng(42)
>>> t = rand.uniform(0, 10, 500)
>>> y = np.ones_like(t)
>>> y[np.abs((t + 1.0)%2.0-1)<0.08] = 1.0 - 0.1
>>> y += 0.01 * rand.standard_normal(len(t))

Compute the transit periodogram on a heuristically determined period grid and find the period with maximum power:

>>> model = BoxLeastSquares(t, y)
>>> results = model.autopower(0.16)
>>> results.period[np.argmax(results.power)]  
2.000412388152837

Compute the periodogram on a user-specified period grid:

>>> periods = np.linspace(1.9, 2.1, 5)
>>> results = model.power(periods, 0.16)
>>> results.power  
array([0.01723948, 0.0643028 , 0.1338783 , 0.09428816, 0.03577543])

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.day
>>> y = y * u.dimensionless_unscaled
>>> model = BoxLeastSquares(t, y)
>>> results = model.autopower(0.16 * u.day)
>>> results.period.unit
Unit("d")
>>> results.power.unit
Unit(dimensionless)

Methods Summary

autoperiod(duration[, minimum_period, ...])

Determine a suitable grid of periods.

autopower(duration[, objective, method, ...])

Compute the periodogram at set of heuristically determined periods.

compute_stats(period, duration, transit_time)

Compute descriptive statistics for a given transit model.

from_timeseries(timeseries[, ...])

Initialize a periodogram from a time series object.

model(t_model, period, duration, transit_time)

Compute the transit model at the given period, duration, and phase.

power(period, duration[, objective, method, ...])

Compute the periodogram for a set of periods.

transit_mask(t, period, duration, transit_time)

Compute which data points are in transit for a given parameter set.

Methods Documentation

autoperiod(duration, minimum_period=None, maximum_period=None, minimum_n_transit=3, frequency_factor=1.0)[source]#

Determine a suitable grid of periods.

This method uses a set of heuristics to select a conservative period grid that is uniform in frequency. This grid might be too fine for some user’s needs depending on the precision requirements or the sampling of the data. The grid can be made coarser by increasing frequency_factor.

Parameters:
durationfloat, array_like, or Quantity [:ref: ‘time’]

The set of durations that will be considered.

minimum_period, maximum_periodfloat or Quantity [:ref: ‘time’], optional

The minimum/maximum periods to search. If not provided, these will be computed as described in the notes below.

minimum_n_transitint, optional

If maximum_period is not provided, this is used to compute the maximum period to search by asserting that any systems with at least minimum_n_transits will be within the range of searched periods. Note that this is not the same as requiring that minimum_n_transits be required for detection. The default value is 3.

frequency_factorfloat, optional

A factor to control the frequency spacing as described in the notes below. The default value is 1.0.

Returns:
periodarray_like or Quantity [:ref: ‘time’]

The set of periods computed using these heuristics with the same units as t.

Notes

The default minimum period is chosen to be twice the maximum duration because there won’t be much sensitivity to periods shorter than that.

The default maximum period is computed as

maximum_period = (max(t) - min(t)) / minimum_n_transits

ensuring that any systems with at least minimum_n_transits are within the range of searched periods.

The frequency spacing is given by

df = frequency_factor * min(duration) / (max(t) - min(t))**2

so the grid can be made finer by decreasing frequency_factor or coarser by increasing frequency_factor.

autopower(duration, objective=None, method=None, oversample=10, minimum_n_transit=3, minimum_period=None, maximum_period=None, frequency_factor=1.0)[source]#

Compute the periodogram at set of heuristically determined periods.

This method calls BoxLeastSquares.autoperiod() to determine the period grid and then BoxLeastSquares.power() to compute the periodogram. See those methods for documentation of the arguments.

compute_stats(period, duration, transit_time)[source]#

Compute descriptive statistics for a given transit model.

These statistics are commonly used for vetting of transit candidates.

Parameters:
periodfloat or Quantity [:ref: ‘time’]

The period of the transits.

durationfloat or Quantity [:ref: ‘time’]

The duration of the transit.

transit_timefloat or Quantity or Time

The mid-transit time of a reference transit.

Returns:
statsdict

A dictionary containing several descriptive statistics:

  • depth: The depth and uncertainty (as a tuple with two

    values) on the depth for the fiducial model.

  • depth_odd: The depth and uncertainty on the depth for a

    model where the period is twice the fiducial period.

  • depth_even: The depth and uncertainty on the depth for a

    model where the period is twice the fiducial period and the phase is offset by one orbital period.

  • depth_half: The depth and uncertainty for a model with a

    period of half the fiducial period.

  • depth_phased: The depth and uncertainty for a model with the

    fiducial period and the phase offset by half a period.

  • harmonic_amplitude: The amplitude of the best fit sinusoidal

    model.

  • harmonic_delta_log_likelihood: The difference in log

    likelihood between a sinusoidal model and the transit model. If harmonic_delta_log_likelihood is greater than zero, the sinusoidal model is preferred.

  • transit_times: The mid-transit time for each transit in the

    baseline.

  • per_transit_count: An array with a count of the number of

    data points in each unique transit included in the baseline.

  • per_transit_log_likelihood: An array with the value of the

    log likelihood for each unique transit included in the baseline.

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 initializing 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(t_model, period, duration, transit_time)[source]#

Compute the transit model at the given period, duration, and phase.

Parameters:
t_modelarray_like, Quantity, or Time

Times at which to compute the model.

periodfloat or Quantity [:ref: ‘time’]

The period of the transits.

durationfloat or Quantity [:ref: ‘time’]

The duration of the transit.

transit_timefloat or Quantity or Time

The mid-transit time of a reference transit.

Returns:
y_modelarray_like or Quantity

The model evaluated at the times t_model with units of y.

power(period, duration, objective=None, method=None, oversample=10)[source]#

Compute the periodogram for a set of periods.

Parameters:
periodarray_like or Quantity [:ref: ‘time’]

The periods where the power should be computed

durationfloat, array_like, or Quantity [:ref: ‘time’]

The set of durations to test

objective{‘likelihood’, ‘snr’}, optional

The scalar that should be optimized to find the best fit phase, duration, and depth. This can be either 'likelihood' (default) to optimize the log-likelihood of the model, or 'snr' to optimize the signal-to-noise with which the transit depth is measured.

method{‘fast’, ‘slow’}, optional

The computational method used to compute the periodogram. This is mainly included for the purposes of testing and most users will want to use the optimized 'fast' method (default) that is implemented in Cython. 'slow' is a brute-force method that is used to test the results of the 'fast' method.

oversampleint, optional

The number of bins per duration that should be used. This sets the time resolution of the phase fit with larger values of oversample yielding a finer grid and higher computational cost.

Returns:
resultsBoxLeastSquaresResults

The periodogram results as a BoxLeastSquaresResults object.

Raises:
ValueError

If oversample is not an integer greater than 0 or if objective or method are not valid.

transit_mask(t, period, duration, transit_time)[source]#

Compute which data points are in transit for a given parameter set.

Parameters:
tarray_like or Quantity [:ref: ‘time’]

Times where the mask should be evaluated.

periodfloat or Quantity [:ref: ‘time’]

The period of the transits.

durationfloat or Quantity [:ref: ‘time’]

The duration of the transit.

transit_timefloat or Quantity or Time

The mid-transit time of a reference transit.

Returns:
transit_maskarray_like

A boolean array where True indicates and in transit point and False indicates and out-of-transit point.