Common Features (starclass.features)

This module contains the features extracted for all targets

Light curve

The primary feature which all classifiers will have available is the light curve of the target in question. This is available in the features dictionary under the key 'lightcurve' and will contain a lightkurve.TessLightCurve object.

Pre-calculated diagnostics

The following features are also provided to the classifiers in the features dictionary- Some of these could be useful in, some may not, but since they are already computed by the previous steps in the TASOC Pipeline, they can be provided to the classifiers “for free”:

  • TESS Magnitude ('tmag').

  • Variance ('variance').

  • RMS on 1 hour timescale ('rms_hour').

  • Point-to-point scatter ('ptp').

Power spectrum (starclass.features.powerspectrum)

Code author: Kristine Kousholt Mikkelsen <201505068@post.au.dk>

Code author: Rasmus Handberg <rasmush@phys.au.dk>

class starclass.features.powerspectrum.powerspectrum(lightcurve, fit_mean=False)[source]

Bases: object

nyquist

Nyquist frequency in Hz.

Type:

float

df

Fundamental frequency spacing in Hz.

Type:

float

standard

Frequency in microHz and power density spectrum sampled from 0 to nyquist with a spacing of df.

Type:

tuple

ls
Type:

astropy.timeseries.LombScargle

Code author: Kristine Kousholt Mikkelsen <201505068@post.au.dk>

Code author: Rasmus Handberg <rasmush@phys.au.dk>

__init__(lightcurve, fit_mean=False)[source]
Parameters:
  • lightcurve (lightkurve.LightCurve) – Lightcurve to estimate power spectrum for.

  • fit_mean (boolean, optional)

alpha_beta(freq)[source]

omega = freq*2*np.pi*1e-6

#w = np.ones_like(self.ls.t) # self.ls.y_err**-2

# Calcultae sums: sx = np.sin(omega*self.ls.t) cx = np.cos(omega*self.ls.t) s = np.sum(self.ls.y * sx) c = np.sum(self.ls.y * cx) cs = np.sum(sx * cx) cc = np.sum(cx * cx)

# Calculate ss on basis of cc and the sum # of the weights, which was calculated in # ImportData: sumWeights = len(self.ls.t) # np.sum(w) ss = sumWeights - cc

# Calculate amplitude and phase: D = ss*cc - cs*cs alpha = (s * cc - c * cs)/D beta = (c * ss - s * cs)/D

copy()[source]

Create copy of power spectrum.

false_alarm_probability(freq)[source]

Calculate Lomb-Scargle false alarm probability for given frequency.

Parameters:

freq (ndarray) – Frequency in microHz.

Returns:

False alarm probability (p-value).

Return type:

ndarray

fundamental_spacing_integral()[source]

Estimate fundamental spacing using the integral of the spectral window function.

fundamental_spacing_minimum()[source]

Estimate fundamental spacing using the first minimum spectral window function.

model(a, b, freq)[source]
optimize_peak(fmax)[source]

Optimize frequency to nearest peak.

Parameters:

fmax (float) – Frequency in microHz.

Returns:

Optimized frequency in microHz.

Return type:

float

plot(ax=None, xlabel='Frequency (muHz)', ylabel=None, style='powerspectrum')[source]
powerspectrum(freq=None, oversampling=1, nyquist_factor=1, scale='power')[source]

Calculate power spectrum for time series.

Parameters:
  • freq (ndarray, optional) – Frequencies to calculate power spectrum for. If set to None, the full frequency range from 0 to nyquist``*``nyquist_factor is calculated.

  • oversampling (float, optional) – Oversampling factor. Default=1.

  • nyquist_factor (float, optional) – Nyquist factor. Default=1.

  • scale (str, optional) – ‘power’, ‘powerdensity’ and ‘amplitude’. Default=’power’.

Returns:

Tuple of two ndarray with frequencies in microHz and corresponding

power in units depending on the scale keyword.

Return type:

tuple

replace_lightcurve(lightcurve)[source]
windowfunction(width=None, oversampling=10)[source]

Spectral window function.

Parameters:
  • width (float, optional) – The width in Hz on either side of zero to calculate spectral window.

  • oversampling (float, optional) – Oversampling factor. Default=10.

Frequency Extraction (starclass.features.freqextr)

Extract frequencies from timeseries.

Code author: Kristine Kousholt Mikkelsen <201505068@post.au.dk>

Code author: Rasmus Handberg <rasmush@phys.au.dk>

starclass.features.freqextr.freqextr(lightcurve, n_peaks=6, n_harmonics=0, hifac=1, ofac=4, snrlim=None, snr_width=None, faplim=0.0027000000000000357, devlim=0.5, conseclim=10, harmonics_list=None, Noptimize=10, optim_max_diff=10, initps=None)[source]

Extract frequencies from timeseries.

The program will perform iterative sine-wave fitting (CLEAN or pre-whitening) using a sum of harmonic functions of the following form:

\[\sum_{i=1}^{N_\mathrm{peaks}} A_i \sin(2\pi\nu_i t + \delta_i) = \sum_{i=1}^{N_\mathrm{peaks}} \alpha_i\sin(2\pi\nu_i t) + \beta_i\cos(2\pi\nu_i t) \, ,\]

where \(\nu_i\), \(A_i\) and \(\delta_i\) denoted the frequency, amplitude and phase of the oscillation.

If n_harmonic is greater than zero, the routine will additionally for each extracted peak extract peaks at the given number of harmonics for each peak. Default is to \(2\nu_i\), \(3\nu_i\) etc., but this can be controlled by the harmonics_list input.

At each iteration, an optimization loop is entered which will go back and re-optimize the previously found peaks in an attempt at minimizing influences between close frequencies. The details of this optimization can be controlled by the parameters Noptimize and optim_max_diff.

Parameters:
  • lightcurve (lightkurve.LightCurve) – Lightcurve to extract frequencies for.

  • n_peaks (int, optional) – Number of frequencies to extract.

  • n_harmonics (int, optional) – Number of harmonics to extract for each frequency.

  • hifac (float, optional) – Nyquist factor.

  • ofac (int, optional) – Oversampling factor used for initial search for peaks in power spectrum.

  • snrlim (float, optional) – Limit on local signal-to-noise ratio above which peaks are considered significant. If set to None no limit is enforced. Default is to not enforce a limit.

  • snr_width (float, optional) – Width in uHz around each peak to estimate signal-to-noise from. Default is 15 frequency steps on either side of the peak.

  • faplim (float, optional) – False alarm probability limit. Peaks with a f.a.p. below this limit are considerd significant. If set to None no limit is enforced. Default is 1-0.9973=0.0027.

  • devlim (float, optional) – If set to None no limit is enforced. Default is 50%.

  • conseclim (int, optional) – Stop after this number of consecutive failed peaks. Default is 10.

  • Noptimize (int, optional) – At each iteration re-optimize. If put to -1, all peaks will be optimized at each iteration. Default is 10.

  • optim_max_diff (float, optional) – Maximal difference in uHz between frequencies to be optimized. Any frequencies futher away than this value from the extracted peak will not be optimized in that iteration. If set to None no limit is enforced. Default is 10 uHz. Please note that this does not take the spectral windowfunction into account, so this value may have to be increased in cases where the windowfunction has significant side-lobes.

  • initps (powerspectrum, optional) – Initial powerspectrum. Should be a powerspectrum calculated from the provided lightcurve. This can be provided if the powerspectrum has already been calculated. If not provided, it is calculated from the provided lightcurve.

Returns:

Table of extracted oscillations.

Return type:

astropy.table.Table

Note

If the height of the peak of one of the harmonics are close to being insignificant, the harmonic may not be found as an harmonic, but will be found later as a peak in it self.

Code author: Kristine Kousholt Mikkelsen <201505068@post.au.dk>

Code author: Rasmus Handberg <rasmush@phys.au.dk>

starclass.features.freqextr.freqextr_table_from_dict(feat, n_peaks=None, n_harmonics=None, flux_unit=None)[source]

Reconstruct freqextr table from features dict.

Please note that not all information can be recreated from features dictionaries. The deviation column will be all NaN, since it can not be recreated, and only n_peaks and n_harmonics will be available in the meta-information.

Parameters:
  • feat (dict) – Dictionary of features.

  • n_peaks (int) – If not provided, it will be determined from feat.

  • n_harmonics (int) – If not provided, it will be determined from feat.

  • flux_unit (astropy.units.Unit) – Unit to put on the amplitude, alpha and beta columns.

Return type:

astropy.table.Table

Code author: Rasmus Handberg <rasmush@phys.au.dk>

starclass.features.freqextr.freqextr_table_to_dict(tab)[source]

Convert freqextr output table to features dictionary.

This will return a dictionary with freq, amp, and phase keys. Please not that this operation does not conserve all information from the table.

Parameters:

tab (astropy.table.Table) – Table extracted by freqextr().

Returns:

Dictionary with frequencies, amplitudes and phases.

Return type:

dict

Code author: Rasmus Handberg <rasmush@phys.au.dk>

starclass.features.freqextr.model(x, a, b, freq)[source]

FliPer (starclass.features.fliper)

This code is the property of L. Bugnet (please see and cite Bugnet et al.,2018).

The user should use the FliPer method to calculate FliPer values from 0.2 ,0.7, 7, 20 and 50 muHz.

Code author: Lisa Bugnet <lisa.bugnet@cea.fr>

Code author: Rasmus Handberg <rasmush@phys.au.dk>

starclass.features.fliper.FliPer(psd)[source]

Compute FliPer values from 0.7, 7, 20, & 50 muHz

Parameters:

psd (powerspectrum object) – Power spectrum of which to calculate the FliPer metrics.

Returns:

Features from FliPer method.

Return type:

dict