wwz_psd (pyleoclim.utils.spectral.wwz_psd)

pyleoclim.utils.spectral.wwz_psd(ys, ts, freq=None, freq_method='log', freq_kwargs=None, tau=None, c=0.001, nproc=8, detrend=False, sg_kwargs=None, gaussianize=False, standardize=False, Neff=3, anti_alias=False, avgs=2, method='Kirchner_numba')[source]

Return the psd of a timeseries using wwz method.

Parameters
  • ys (array) – a time series, NaNs will be deleted automatically

  • ts (array) – the time points, if ys contains any NaNs, some of the time points will be deleted accordingly

  • freq (array) – vector of frequency

  • freq_method (str, {'log', 'lomb_scargle', 'welch', 'scale', 'nfft'}) –

    Method to generate the frequency vector if not set directly. The following options are avialable:

    • ’log’ (default)

    • ’lomb_scargle’

    • ’welch’

    • ’scale’

    • ’nfft’

    See pyleoclim.utils.wavelet.make_freq_vector() for details

  • freq_kwargs (dict) – Arguments for the method chosen in freq_method. See specific functions in pyleoclim.utils.wavelet for details

  • tau (array) – the evenly-spaced time vector for the analysis, namely the time shift for wavelet analysis

  • c (float) – the decay constant that will determine the analytical resolution of frequency for analysis, the smaller the higher resolution; the default value 1e-3 is good for most of the spectral analysis cases

  • nproc (int) – the number of processes for multiprocessing

  • detrend (str, {None, 'linear', 'constant', 'savitzy-golay'}) –

    available methods for detrending, including

    • None: the original time series is assumed to have no trend;

    • ’linear’: a linear least-squares fit to ys is subtracted;

    • ’constant’: the mean of ys is subtracted

    • ’savitzy-golay’: ys is filtered using the Savitzky-Golay filters and the resulting filtered series is subtracted from y.

  • sg_kwargs (dict) – The parameters for the Savitzky-Golay filters. See pyleoclim.utils.filter.savitzky_golay() for details.

  • gaussianize (bool) – If True, gaussianizes the timeseries

  • standardize (bool) – If True, standardizes the timeseries

  • method (string, {'Foster', 'Kirchner', 'Kirchner_f2py', 'Kirchner_numba'}) –

    available specific implementation of WWZ, including

    • ’Foster’: the original WWZ method;

    • ’Kirchner’: the method Kirchner adapted from Foster;

    • ’Kirchner_f2py’: the method Kirchner adapted from Foster, implemented with f2py for acceleration;

    • ’Kirchner_numba’: the method Kirchner adapted from Foster, implemented with Numba for acceleration (default);

  • Neff (int) – effective number of points

  • anti_alias (bool) – If True, uses anti-aliasing

  • avgs (int) – flag for whether spectrum is derived from instantaneous point measurements (avgs<>1) OR from measurements averaged over each sampling interval (avgs==1)

Returns

res – a namedtuple that includes below items

psdarray

power spectral density

freqarray

vector of frequency

Return type

namedtuple

See also

pyleoclim.utils.spectral.periodogram

Estimate power spectral density using a periodogram

pyleoclim.utils.spectral.mtm

Retuns spectral density using a multi-taper method

pyleoclim.utils.spectral.lomb_scargle

Return the computed periodogram using lomb-scargle algorithm

pyleoclim.utils.spectral.welch

Estimate power spectral density using the Welch method

pyleoclim.utils.filter.savitzy_golay

Filtering using Savitzy-Golay

pyleoclim.utils.tsutils.detrend

Detrending method

References

  • Foster, G. (1996). Wavelets for period analysis of unevenly sampled time series. The Astronomical Journal, 112(4), 1709-1729.

  • Kirchner, J. W. (2005). Aliasin in 1/f(alpha) noise spectra: origins, consequences, and remedies. Physical Review E covering statistical, nonlinear, biological, and soft matter physics, 71, 66110.