Periodograms - periodogram estimation
#
DSP.Periodograms.arraysplit — Function
arraysplit(s, n, m)
Split an array into arrays of length n with overlapping regions of length m. Iterating or indexing the returned AbstractVector always yields the same Vector with different contents.
#
DSP.Periodograms.periodogram — Method
periodogram(s; onesided=eltype(s)<:Real, nfft=nextfastfft(n), fs=1, window=nothing)
Computes periodogram of a signal by FFT and returns a Periodogram object.
For real signals, the two-sided periodogram is symmetric and this function returns a one-sided (real only) periodogram by default. A two-sided periodogram can be obtained by setting onesided=false.
nfft specifies the number of points to use for the Fourier transform. If length(s) < nfft, then the input is padded with zeros. By default, nfft is the closest size for which the Fourier transform can be computed with maximal efficiency.
fs is the sample rate of the original signal, and window is an optional window function or vector to be applied to the original signal before computing the Fourier transform. The computed periodogram is normalized so that the area under the periodogram is equal to the uncentered variance (or average power) of the original signal.
#
DSP.Periodograms.welch_pgram — Function
welch_pgram(s, n=div(length(s), 8), noverlap=div(n, 2); onesided=eltype(s)<:Real, nfft=nextfastfft(n), fs=1, window=nothing)
Computes the Welch periodogram of a signal s based on segments with n samples with overlap of noverlap samples, and returns a Periodogram object. For a Bartlett periodogram, set noverlap=0. See periodogram for description of optional keyword arguments.
#
DSP.Periodograms.spectrogram — Function
spectrogram(s, n=div(length(s), 8), noverlap=div(n, 2); onesided=eltype(s)<:Real, nfft=nextfastfft(n), fs=1, window=nothing)
Computes the spectrogram of a signal s based on segments with n samples with overlap of noverlap samples, and returns a Spectrogram object. See periodogram for description of optional keyword arguments.
#
DSP.Periodograms.stft — Function
stft(s, n=div(length(s), 8), noverlap=div(n, 2); onesided=eltype(s)<:Real, nfft=nextfastfft(n), fs=1, window=nothing)
Computes the STFT of a signal s based on segments with n samples with overlap of noverlap samples, and returns a matrix containing the STFT coefficients. See periodogram for description of optional keyword arguments.
#
DSP.Periodograms.periodogram — Method
periodogram(s::AbstractMatrix; nfft=nextfastfft(size(s)), fs=1, radialsum=false, radialavg=false)
Computes periodogram of a 2-d signal by FFT and returns a Periodogram2 object.
Returns a 2-d periodogram by default. A radially summed or averaged periodogram is returned as a Periodogram object if radialsum or radialavg is true, respectively.
nfft specifies the number of points to use for the Fourier transform. If size(s) < nfft, then the input is padded with zeros. By default, nfft is the closest size for which the Fourier transform can be computed with maximal efficiency. fs is the sample rate of the original signal in both directions.
For radialsum=true the value of power[k] is proportional to |^2]. For radialavg=true it is proportional to |^2]. The computation of \|k'\| takes into account non-square signals by scaling the coordinates of the wavevector accordingly.
#
DSP.Periodograms.freq — Function
freq(p)
Returns the frequency bin centers for a given Periodogram, Spectrogram, CrossPowerSpectra, or Coherence object.
Returns a tuple of frequency bin centers for a given Periodogram2 object.
#
DSP.Periodograms.power — Function
power(p)
For a Periodogram, returns the computed power at each frequency as a Vector.
For a Spectrogram, returns the computed power at each frequency and time bin as a Matrix. Dimensions are frequency × time.
For a CrossPowerSpectra, returns the pairwise power between each pair of channels at each frequency. Dimensions are channel x channel x frequency.
#
Base.Libc.time — Function
time(p)
Returns the time bin centers for a given Spectrogram object.
#
DSP.Periodograms.coherence — Function
coherence(c::Coherence)
Given an Coherence object, returns an n_channels x n_channels x length(freq(c)) array consisting of the pairwise coherences between each channel for each frequency.
Multitaper periodogram estimation
#
DSP.Periodograms.mt_pgram — Function
mt_pgram(s; onesided=eltype(s)<:Real, nfft=nextfastfft(n), fs=1, nw=4, ntapers=iceil(2nw)-1, window=dpss(length(s), nw, ntapers))
mt_pgram(signal::AbstractVector, config::MTConfig)
Computes the multitaper periodogram of a signal s.
If window is not specified, the signal is tapered with ntapers discrete prolate spheroidal sequences with time-bandwidth product nw. Each sequence is equally weighted; adaptive multitaper is not (yet) supported.
If window is specified, each column is applied as a taper. The sum of periodograms is normalized by the total sum of squares of window.
Returns a Periodogram.
#
DSP.Periodograms.mt_pgram! — Function
mt_pgram!(output, s::AbstractVector{T}; onesided::Bool=eltype(s)<:Real,
nfft::Int=nextfastfft(length(s)), fs::Real=1,
nw::Real=4, ntapers::Int=ceil(Int, 2nw)-1,
window::Union{AbstractMatrix,Nothing}=nothing) where T<:Number
mt_pgram!(output::AbstractVector, signal::AbstractVector, config::MTConfig) -> Periodogram
Computes a multitapered periodogram, storing the output in output. Arguments:
-
signal::AbstractVector: should be of lengthconfig.n_samples -
output::AbstractVector: should be of lengthlength(config.freq)
Optionally pass an MTConfig object to preallocate temporary variables and choose configuration settings; otherwise, keyword arguments may be passed to choose those settings.
Returns a Periodogram.
#
DSP.Periodograms.mt_spectrogram — Function
mt_spectrogram(signal::AbstractVector{T}, n::Int=length(s) >> 3,
n_overlap::Int=n >> 1; fs=1,
onesided::Bool=T <: Real, kwargs...) where {T}
mt_spectrogram(signal::AbstractVector, config::MTSpectrogramConfig)
Compute a multitaper spectrogram, returning a Spectrogram object. Optionally pass a MTSpectrogramConfig object; otherwise, any additional keyword arguments accepted by MTConfig may be passed to configure the tapering.
Returns a Spectrogram.
See also mt_spectrogram!.
#
DSP.Periodograms.mt_spectrogram! — Function
mt_spectrogram!(output, signal::AbstractVector{T}, n::Int=length(signal) >> 3,
n_overlap::Int=n >> 1; fs=1, onesided::Bool=T <: Real, kwargs...) where {T}
mt_spectrogram!(destination::AbstractMatrix, signal::AbstractVector, config::MTSpectrogramConfig)
Computes a multitaper spectrogram using the parameters specified in config. Arguments:
-
destination:length(config.mt_config.freq)xlength(config.time)matrix. This can be created byDSP.allocate_output(config). -
signal: vector of lengthconfig.n_samples -
config: optionally, pass anMTSpectrogramConfigobject to hold temporary variables and configuration settings. Otherwise, settings arguments may be passed directly.
Returns a Spectrogram.
See also mt_spectrogram.
#
DSP.Periodograms.mt_cross_power_spectra — Function
mt_cross_power_spectra(signal::AbstractMatrix{T}; fs=1, kwargs...) where {T}
mt_cross_power_spectra(signal::AbstractMatrix, config::MTCrossSpectraConfig)
Computes multitapered cross power spectra between channels of a signal. Arguments:
-
signal:n_channelsxn_samples -
Optionally pass an
MTCrossSpectraConfigobject to preallocate temporary variables
and choose configuration settings. Otherwise, any keyword arguments accepted by MTCrossSpectraConfig may be passed here.
Produces a CrossPowerSpectra object holding the n_channels x n_channels x n_frequencies output array (accessed by power) and the corresponding frequencies (accessed by freq).
See also mt_cross_power_spectra! and MTCrossSpectraConfig.
#
DSP.Periodograms.mt_cross_power_spectra! — Function
mt_cross_power_spectra!(output, signal::AbstractMatrix; fs=1, kwargs...)
mt_cross_power_spectra!(output, signal::AbstractMatrix, config::MTCrossSpectraConfig)
Computes multitapered cross power spectra between channels of a signal. Arguments:
-
output:n_channelsxn_channelsxlength(config.freq). Can be created byDSP.allocate_output(config). -
signal:n_channelsxn_samples -
config:MTCrossSpectraConfig{T}: optionally pass aMTCrossSpectraConfigto preallocate temporary and choose configuration settings. Otherwise, one may pass any keyword arguments accepted by this object.
Produces a CrossPowerSpectra object holding the n_channels x n_channels x n_frequencies output array and the corresponding frequencies (accessed by freq).
See also mt_cross_power_spectra and MTCrossSpectraConfig.
#
DSP.Periodograms.mt_coherence — Function
mt_coherence(signal::AbstractMatrix{T}; fs=1, freq_range = nothing, demean=false, kwargs...) where T
mt_coherence(signal::AbstractMatrix, config::MTCoherenceConfig)
Arguments:
-
signal:n_channelsxn_samplesmatrix -
Optionally pass an
MTCoherenceConfigto pre-allocate temporary variables and choose configuration settings, otherwise, seeMTCrossSpectraConfigfor the meaning of the keyword arguments.
Returns a Coherence object.
See also mt_coherence and MTCoherenceConfig.
#
DSP.Periodograms.mt_coherence! — Function
mt_coherence!(output, signal::AbstractMatrix; fs=1, freq_range=nothing, demean=false, kwargs...)
mt_coherence!(output, signal::AbstractMatrix, config::MTCoherenceConfig)
Computes the pairwise coherences between channels.
-
output:n_channelsxn_channelsmatrix -
signal:n_samplesxn_channelsmatrix -
config: optional configuration object that pre-allocates temporary variables and choose settings.
Returns a Coherence object.
See also mt_coherence and MTCoherenceConfig.
Configuration objects
#
DSP.Periodograms.MTConfig — Type
MTConfig{T}(n_samples; fs=1,
nfft = nextpow(2, n_samples),
window = nothing,
nw = 4,
ntapers = 2 * nw - 1,
taper_weights = fill(1/ntapers, ntapers),
onesided::Bool=T<:Real,
fft_flags = FFTW.MEASURE)
Creates a config object which holds the configuration state and temporary variables used in multitaper computations, e.g. mt_pgram!, mt_spectrogram, MTSpectrogramConfig, MTCrossSpectraConfig, and MTCoherenceConfig.
An MTConfig can be re-used between computations as long as none of the input arguments change.
-
n_samples: the number of samples to be used as input when computing multitaper periodograms with this configuration. Used for pre-allocating temporary buffers. -
fs: the number of samples per second of the input signal -
nfft: length of input vector to the FFT; ifnfft > n_samples, then the input signal will be zero-padded until it is of lengthnfft. -
window: window function to use for tapering. If left at the default ofnothing,windowwill be set todpss(n_samples, nw, ntapers). -
ntapers: the number of tapers to use. -
taper_weights = fill(1/ntapers, ntapers): how to weight the contribution of each taper. The default setting is to simply average them. -
onesided: whether or not to compute a "one-sided" FFT by using that real signal data yields conjugate-symmetry in Fourier space. -
fft_flags: flags to control how the FFT plan is generated.
#
DSP.Periodograms.MTSpectrogramConfig — Type
MTSpectrogramConfig(n_samples, mt_config::MTConfig{T}, n_overlap_samples) where {T}
MTSpectrogramConfig{T}(n_samples, samples_per_window, n_overlap_samples; fs=1, kwargs...) where {T}
Creates a MTSpectrogramConfig which holds configuration and temporary variables for mt_spectrogram and mt_spectrogram!. Any keyword arguments accepted by MTConfig may be passed here, or an MTConfig object itself.
#
DSP.Periodograms.MTCrossSpectraConfig — Type
MTCrossSpectraConfig{T}(n_channels, n_samples; fs=1, demean=false, freq_range=nothing,
ensure_aligned = T == Float32 || T == Complex{Float32}, kwargs...) where {T}
MTCrossSpectraConfig(n_channels, mt_config::MTConfig{T}; demean=false, freq_range=nothing,
ensure_aligned = T == Float32 || T == Complex{Float32})
Creates a configuration object used for mt_cross_power_spectra and mt_cross_power_spectra!.
-
n_channels: the number of channels of the input. -
n_samples: the number of samples for each channel of the input. -
demean: iftrue, the channelwise mean will be subtracted from the input signals before the cross spectral powers are computed. -
freq_range: ifnothing, all frequencies are retained. Otherwise, only frequencies betweenfirst(freq_range)andlast(freq_range)are retained. -
ensure_aligned = T == Float32 || T == Complex{Float32}: perform an extra copy to ensure that the FFT output is memory-aligned -
Additionally, either pass an
MTConfigobject, or keyword arguments such asfsaccepted byMTConfig.
Returns a CrossPowerSpectra object.
#
DSP.Periodograms.MTCoherenceConfig — Type
MTCoherenceConfig{T}(n_channels, n_samples; fs=1, demean=false, freq_range=nothing, kwargs...) where T
MTCoherenceConfig(cs_config::MTCrossSpectraConfig{T}) where {T}
MTCoherenceConfig(n_channels, mt_config::MTConfig{T}; demean=false, freq_range=nothing,
ensure_aligned = T == Float32 || T == Complex{Float32}) where {T}
Creates a configuration object for coherences from a MTCrossSpectraConfig. Provides a helper method with the same arugments as MTCrossSpectraConfig to construct the MTCrossSpectraConfig object.
See also mt_coherence and mt_coherence!.