Function index

Public types and functions

Structs

Multitaper.MTSpectrumType

The multitaper spectrum is given as a MTSpectrum structure which holds

  • frequency (f),
  • spectrum (S),
  • phase (optional),
  • chosen values of the multitaper time bandwidth product etc of type MTParameters (params)
  • eigencoefficients (coef, optional),
  • Ftest values (Fpval, optional),
  • jackknife output (jkvar, optional), and
  • Tsquared test results (Tsq_pval, optional).
Multitaper.MTCoherenceType

The multitaper coherence structure, MTCoherence, holds

  • frequency (f),
  • coherence (coh),
  • phase (phase),
  • eigencoefficients (coef, optional),
  • jackknife output (jkvar, optional), and
  • Tsquared test results (Tsq, optional).
Multitaper.MTTransferFunctionType

The multitaper transfer function is given as a MTTransferFunction structure which holds

  • frequency (f),
  • transfer function (transf),
  • phase (phase),
  • MTParameters (params),
  • eigencoefficients (EigenCoefficientfs, optional),
  • jackknife output (jkvar, optional), and
  • Tsquared test results (Tsq, optional).
Multitaper.EigenCoefficientType

The EigenCoefficient structure holds

  • multitaper eigencoefficients (coef) and, optionally,
  • adaptive weights (wts)
Multitaper.MTParametersType

Multitaper parameters MTParameters struct contains

  • time bandwidth (NW) as Float,
  • number of tapers (K),
  • number of samples (N),
  • sampling rate (dt) in temporal units (e.g. seconds),
  • padded length (M),
  • number of segments (nsegments), and
  • overlap if nsegments is greater than 1, nothing otherwise.
Multitaper.DemodulateType

Multitaper complex demodulates are held in the Demodulate struct which contains

  • magnitude (Mag),
  • phase (Phase)
Multitaper.MTAutocorrelationFunctionType

The multitaper autocorrelation function is given in the MTAutocorrelationFunction structure which holds

  • lags (lags),
  • autocorrelation function (acf),
  • MTParameters (params)
Multitaper.MTAutocovarianceFunctionType

The multitaper autocovariance function is given in the MTAutocovarianceFunction structure, which holds

  • lags (lags),
  • autocovariance function (acvf),
  • MTParameters (params)
Multitaper.MTCepstrumType

The multitaper cepstrum is given in the MTCepstrum structure which holds

  • lags (lags),
  • cepstrum (ceps),
  • MTParameters (params)
Multitaper.MtCrossCovarianceFunctionType

The multitaper cross-covariance function is given in the MtCrossCovarianceFunction structure which holds

  • lags (lags),
  • cross-covariance function (acvf),
  • MTParameters (params)
Multitaper.MTCrossCorrelationFunctionType

The multitaper cross-correlation function is given as a MTCrossCorrelationFunction structure which holds

  • lags (lags),
  • cross-correlation function (acf),
  • MTParameters (params)

Taper functions

Multitaper.dpss_tapersFunction
dpss_tapers(n,w,k,tap_or_egval)

Simply compute discrete prolate spheroidal sequence tapers, eigenvalues

...

Arguments

  • n::Int64: Length of the taper

  • nw::Float64: Time-bandwidth product

  • k::Int64: Number of tapers

  • tap_or_egval::Symbol = :tap: Either :tap, :egval, or :both

...

...

Outputs

  • vv::Vector{Float64}: The matrix of eigenvalues, if taporegval is set to :tap

  • dpss_eigval: Struct conaining the dpss tapers

...

Multitaper.mdslepianFunction
mdslepian(w, k, t)

Generalized prolate spheroidal sequences for the 1D missing data problem

...

Arguments

Positional Arguments

  • w::Float64: the bandwidth

  • k::Int64: number of Slepian tapers, must be <=2bwlength(x)

  • t::Vector{Int64}: vector containing the time indices

...

...

Outputs

  • lambda,u::Tuple{Vector{Float64}, Vector{Float64}}: tuple containing the

concentrations and the tapers

...

See also: mdmultispec, gpss

Multitaper.gpssFunction
gpss(w, k, t, f; <keyword arguments>)

Generalized prolate spheroidal sequences on an unequal grid

...

Arguments

Positional Arguments

  • w::Float64: the bandwidth

  • k::Int64: number of Slepian tapers, must be <=2bwlength(x)

  • t::Vector{Int64}: vector containing the time indices

  • f::Float64: frequency at which the tapers are to be computed

Keyword Arguments

  • beta::Float64 = 0.5: analysis half-bandwidth (similar to Nyquist rate)

...

...

Outputs

  • lambda::Vector{Float64} the concentrations of the generalized prolate spheroidal

sequences

  • u::Matrix{Float64} the matrix containing the sequences themselves

  • R the Cholesky factor for the generalized eigenvalue problem

...

This function is currently not exported, use Multitaper.gpss.

See also: mdmultispec, mdslepian

Multitaper Spectrum

Multitaper.multispecFunction
multispec(S1; <keyword arguments>)

Computes univariate multitaper spectra with a handful of extra gadgets.

...

Arguments

  • S1::Vector{T} where T<:Float64: the vector containing the time series
  • NW::Float64 = 4.0: time-bandwidth product of estimate
  • K::Int64 = 6: number of slepian tapers, must be <= 2*NW
  • dt::Float64: sampling rate in time units
  • ctr::Bool: whether or not to remove the mean before computing the multitaper spectrum
  • pad::Float64 = 1.0: factor by which to pad the series, i.e. spectrum length will be pad times length of the time series.
  • dpVec::Union{Matrix{Float64},Nothing} = nothing: Matrix of dpss's, if they have been precomputed
  • egval::Union{Vector{Float64},Nothing} = nothing: Vector of concentratins of said dpss's
  • guts::Bool = false: whether or not to return the eigencoefficients in the output struct
  • a_weight::Bool = true: whether or not to use adaptive weighting
  • Ftest::Bool = true: Compute the F-test p-value
  • highres::Bool = false: Whether to return a "high resolution" spectrum estimate
  • jk::Bool = true: Compute jackknifed confidence intervals
  • Tsq::Union{Vector{Int64},Vector{Vector{Int64}},Nothing} = nothing: which frequency indices to compute the T-squared test for multiple line components. Defaults to none.

...

...

Outputs

  • MTSpectrum struct containing the spectrum

...

See also: dpss_tapers, MTSpectrum, mdmultispec, mdslepian

multispec(S1, S2; <keyword arguments>)

Computes multitaper cross-spectrum or coherence when given two time series with same sampling.

...

Arguments

  • S1::Union{Vector{T},EigenCoefficient} where T<:Number: the vector containing the first time

series

  • S2::Union{Vector{P},EigenCoefficient} where P<:Number: the vector containing the second

time series

  • outp::Symbol: output can be either :coh for coherence, :spec for cross-spectrum,

or :transf for transfer function

  • NW::Float64 = 4.0: time-bandwidth product of estimate

  • K::Int64 = 6: number of slepian tapers, must be <= 2*NW

  • offset::Union{Float64,Int64} = 0 set to nonzero value if offset coherence or

cross-spectrum is desired. If Float64 is used, this will be converted to nearest FFT bin.

  • dt::Float64: sampling rate in time units

  • ctr::Bool: whether or not to remove the mean before computing the multitaper

spectrum

  • pad::Float64 = 1.0: factor by which to pad the series, i.e. spectrum length will

be pad times length of the time series.

  • dpVec::Union{Matrix{Float64},Nothing} = nothing: Matrix of dpss's, if they have

been precomputed

  • guts::Bool = false: whether or not to return the eigencoefficients in the output

struct

  • jk::Bool = true: Compute jackknifed confidence intervals

  • Tsq::Union{Vector{Int64},Vector{Vector{Int64}},Nothing} = nothing: which

frequency indices to compute the T-squared test for multiple line components. Defaults to none.

  • alph::Float64 = 0.05: significance cutoff for the Tsquared test

...

...

Outputs

  • MTSpectrum, MTCoherence, or MTTransferFunction struct containing the spectrum, coherence or

transfer function, depending on the selection of outp input.

...

See also: dpss_tapers, MTSpectrum, mdmultispec, mdslepian

multispec(S1; <keyword arguments>)

Multivariate version of the multispec call, data are in the columns of a matrix ...

Arguments

  • S1::Matrix{T} where T<:Float64: the vector containing the first time series
  • outp::Symbol: output can be either :coh for coherence, :justspeccs to compute just the spectra, or :cross for cross-spectra
  • NW::Float64 = 4.0: time-bandwidth product of estimate
  • K::Int64 = 6: number of slepian tapers, must be <= 2*NW
  • dt::Float64: sampling rate in time units
  • ctr::Bool: whether or not to remove the mean before computing the multitaper spectrum
  • pad::Float64 = 1.0: factor by which to pad the series, i.e. spectrum length will be pad times length of the time series.
  • dpVec::Union{Matrix{Float64},Nothing} = nothing: Matrix of dpss's, if they have been precomputed
  • guts::Bool = false: whether or not to return the eigencoefficients in the output struct
  • a_weight::Bool = true: whether or not to adaptively weight the spectra
  • jk::Bool = false: Compute jackknifed confidence intervals
  • Ftest:Bool = false: Compute F-test for line components
  • Tsq::Union{Vector{Int64},Vector{Vector{Int64}},Nothing} = nothing: which frequency indices to compute the T-squared test for multiple line components. Defaults to none.
  • alph::Float64 = 0.05: significance cutoff for the Tsquared test

...

...

Outputs

  • Tuple{Vector{MTSpectrum},Vector{P},Union{Float64,Vector{Float64}}} where P = Union{MTCoherence,MTSpectrum}

struct containing the spectra, coherence or crossspectra, and Tsquared test p-values. Ouput of middle arg depends on the selection of outp input. ...

See also: dpss_tapers, MTSpectrum, mdmultispec, mdslepian

Multitaper.mdmultispecFunction
mdmultispec(tt, x; <keyword arguments>)

Multitaper power spectrum estimation for time series with missing data (gaps)

...

Arguments

Positional Arguments

  • tt::Vector{T} where T<:Real: the vector containing the time indices

  • x::Vector{P} where P<:Number: data vector

Keyword Arguments

  • bw = 5/length(tt): bandwidth of estimate

  • k::Int64 = 2*bw*length(x)-1: number of Slepian tapers, must be `<=

2bwlength(x)`

  • dt::T = tt[2]-tt[1]: sampling rate in time units

  • nz = 0.0: zero padding factor

  • Ftest::Bool = true: Compute the F-test p-value

  • jk::Bool = true: Compute jackknifed confidence intervals

  • dof::Bool = false: Compute degrees of freedom for the adaptively weighted

spectrum estimate

  • lambdau::Union{Tuple{Array{Float64,1},Array{Float64,2}},Nothing} = nothing:

Slepians, if precomputed ...

...

Outputs

  • pkg::MTSpectrum struct containing the spectrum

  • nu1::Vector{Float64} optional vector containing the degrees of freedom, given

if the dof kwarg is set to true.

...

See also: multispec, mdslepian

mdmultispec(tt, x, y; <keyword arguments>)

Multitaper coherence estimation for time series with missing data (gaps)

...

Arguments

Positional Arguments

  • tt::Vector{T} where T<:Real: the vector containing the time indices

  • x::Vector{P} where P<:Number: data vector 1

  • y::Vector{Q} where Q<:Number: data vector 2

Keyword Arguments

  • bw = 5/length(tt): bandwidth of estimate

  • k::Int64 = 2*bw*length(x)-1: number of Slepian tapers, must be `<=

2bwlength(x)`

  • dt = tt[2]-tt[1]: sampling rate in time units

  • nz = 0.0: zero padding factor

  • Ftest::Bool = true: Compute the F-test p-value

  • jk::Bool = true: Compute jackknifed confidence intervals

  • lambdau::Union{Tuple{Array{Float64,1},Array{Float64,2}},Nothing} = nothing:

Slepians, if precomputed

...

...

Outputs

  • pkg::MTCoherence struct containing the coherence

...

See also: multispec, mdslepian

mdmultispec(tt, x; <keyword arguments>)

Multitaper coherence estimation for multiple time series with the same missing data (gaps)

...

Arguments

Keyword Arguments

  • tt::Vector{T} where T<:Real: the vector containing the time indices

  • x::Matrix{P} where P<:Number: time series in the columns of a matrix

Positional Arguments

  • bw = 5/length(tt): bandwidth of estimate

  • k::Int64 = 2*bw*length(x)-1: number of Slepian tapers, must be `<=

2bwlength(x)`

  • dt = tt[2]-tt[1]: sampling rate in time units

  • nz = 0.0: zero padding factor

  • Ftest::Bool = false: Compute the F-test p-value

  • jk::Bool = true: Compute jackknifed confidence intervals

  • lambdau::Union{Tuple{Array{Float64,1},Array{Float64,2}},Nothing} = nothing:

Slepians, if precomputed

...

...

Outputs

  • Tuple{Vector{MTSpectrum},Matrix{MTCoherence},Nothing} struct containing the spectra,

coherences, and T^2 test significances (currently set to return nothing)

...

See also: multispec, mdslepian

Multitaper.welchFunction
welch(S1, nsegments; <keyword arguments>)

Computes univariate multitaper Welch specrum starting with an input time series.

...

Arguments

  • S1:Vector{T} where T<:Number: the vector containing the time series

  • nsegments::Int64: the number of segments into which to divide the time series

  • overlap::Float64 = 0.5: number between 0 and 1 which accounts for the amount of

overlap between the segments

  • NW::Float64 = 4.0: time-bandwidth product of estimate

  • K::Int64 = 6: number of slepian tapers, must be <= 2*NW

  • dt::Float64: sampling rate in time units

  • ctr::Bool: whether or not to remove the mean before computing the multitaper

spectrum

  • pad::Float64 = 1.0: factor by which to pad the series, i.e. spectrum length will

be pad times length of the time series.

  • dpVec::Union{Matrix{Float64},Nothing} = nothing: Matrix of dpss's, if they have

been precomputed

  • egval::Union{Vector{Float64},Nothing} = nothing: Vector of concentratins of said

dpss's

  • guts::Bool = false: whether or not to return the eigencoefficients in the output

struct

  • a_weight::Bool = true: whether or not to use adaptive weighting

...

...

Outputs

  • MTSpectrum struct, depending on the selection of outp above

  • Float64 containing the effective bandwidth

...

See also: multispec

Time-domain Statistics

Multitaper.mt_acfFunction
mt_acf(S)

Computes univariate multitaper autocorrelation function. Inputs a MTSpectrum struct.

...

Arguments

  • S::MTSpectrum: the vector containing the result of an univariate call to multispec

...

...

Outputs

  • MTAutocorrelationFunction struct containing the autocorrelation function

...

See also: multispec

mt_acf(S1; <keyword arguments>)

Computes univariate multitaper autocorrelation function starting with an input time series. ...

Arguments

  • S1::Vector{T} where T<:Number: the vector containing the time series
  • NW::Float64 = 4.0: time-bandwidth product of estimate
  • K::Int64 = 6: number of slepian tapers, must be <= 2*NW
  • dt::Float64: sampling rate in time units
  • ctr::Bool: whether or not to remove the mean before computing the multitaper spectrum
  • pad::Float64 = 1.0: factor by which to pad the series, i.e. spectrum length will be pad times length of the time series.
  • dpVec::Union{Matrix{Float64},Nothing} = nothing: Matrix of dpss's, if they have been precomputed
  • egval::Union{Vector{Float64},Nothing} = nothing: Vector of concentratins of said dpss's
  • a_weight::Bool = true: whether or not to use adaptive weighting

...

...

Outputs

  • MTAutocorrelationFunction struct containing the autocorrelation function

...

See also: multispec

Multitaper.mt_acvfFunction
mt_acvf(S)

Computes univariate multitaper autocovariance function. Inputs a MTSpectrum struct.

...

Arguments

  • S::MTSpectrum: the vector containing the result of an univariate call to multispec

...

...

Outputs

  • MTAutocovarianceFunction struct containing the autocovariance function.

...

See also: multispec

mt_acvf(S1; <keyword arguments>)

Computes univariate multitaper autocovariance function starting with an input time series. ...

Arguments

  • S1::Vector{T} where T<:Number: the vector containing the time series
  • NW::Float64 = 4.0: time-bandwidth product of estimate
  • K::Int64 = 6: number of slepian tapers, must be <= 2*NW
  • dt::Float64: sampling rate in time units
  • ctr::Bool: whether or not to remove the mean before computing the multitaper spectrum
  • pad::Float64 = 1.0: factor by which to pad the series, i.e. spectrum length will be pad times length of the time series.
  • dpVec::Union{Matrix{Float64},Nothing} = nothing: Matrix of dpss's, if they have been precomputed
  • egval::Union{Vector{Float64},Nothing} = nothing: Vector of concentratins of said dpss's
  • a_weight::Bool = true: whether or not to use adaptive weighting

...

...

Outputs

  • MTAutocovarianceFunction struct containing the autocovariance function

...

See also: multispec

Multitaper.mt_ccvfFunction
mt_ccvf(S; <keyword arguments>)

Computes univariate multitaper cross-covariance/cross-correlation function. Inputs a MTSpectrum struct.

...

Arguments

  • S::MTSpectrum: the vector containing the result of an multiivariate call to multispec
  • typ::Symbol = :ccvf: whether to compute cross-correlation function (:ccf) or cross-covariance function (:ccvf)

...

...

Outputs

  • MtCrossCovarianceFunction, MTCrossCorrelationFunction depending on the selection of typ input above.

...

See also: multispec

mt_ccvf(S1, S2; <keyword arguments>)

Computes bivariate multitaper cross-covariance/cross-correlation function from two time series

...

Arguments

  • S1::Vector{T} where T<:Number: the vector containing the first time series
  • S2::Vector{T} where T<:Number: the vector containing the second time series
  • typ::Symbol: whether to compute cross-covariance function (:ccvf), or cross-correlation function (:ccf)
  • NW::Float64 = 4.0: time-bandwidth product of estimate
  • K::Int64 = 6: number of slepian tapers, must be <= 2*NW
  • dt::Float64: sampling rate in time units
  • ctr::Bool: whether or not to remove the mean before computing the multitaper spectrum
  • pad::Float64 = 1.0: factor by which to pad the series, i.e. spectrum length will be pad times length of the time series.
  • dpVec::Union{Matrix{Float64},Nothing} = nothing: Matrix of dpss's, if they have been precomputed

...

...

Outputs

  • MtCrossCovarianceFunction struct, depending on the selection of typ input above.

...

See also: multispec

Multitaper.mt_cepstrumFunction
mt_cepstrum(S)

Computes multitaper cepstrum. Inputs a MTSpectrum struct.

...

Arguments

  • S::MTSpectrum: the vector containing the result of an univariate call to multispec

...

...

Outputs

  • MTCepstrum struct containing the cepstrum

...

See also: multispec

mt_cepstrum(S1; <keyword arguments>)

Computes multitaper cepstrum starting with an input time series. ...

Arguments

  • S1::Vector{T} where T<:Number: the vector containing the time series
  • NW::Float64 = 4.0: time-bandwidth product of estimate
  • K::Int64 = 6: number of slepian tapers, must be <= 2*NW
  • dt::Float64: sampling rate in time units
  • ctr::Bool: whether or not to remove the mean before computing the multitaper spectrum
  • pad::Float64 = 1.0: factor by which to pad the series, i.e. spectrum length will be pad times length of the time series.
  • dpVec::Union{Matrix{Float64},Nothing} = nothing: Matrix of dpss's, if they have been precomputed
  • egval::Union{Vector{Float64},Nothing} = nothing: Vector of concentratins of said dpss's
  • a_weight::Bool = true: whether or not to use adaptive weighting

...

...

Outputs

  • MTCepstrum struct containing the cepstrum

...

See also: multispec

Multitaper.demodulateFunction
demodulate(x, f0, NW, blockLen; <keyword arguments>)

Multitaper complex demodulation

...

Arguments

  • x::Vector{T} where T<:Float64: the vector containing the time series
  • f0::Float64: center frequency for the complex demodulate
  • NW::Float64: the time bandwidth product
  • blockLen::Int64: the length of the blocks to use
  • wrapphase::Bool = true: whether or not to unwrap the phase
  • dt::Float64 = 1.0: the sampling rate of the series, in time units
  • basetime::Float64 = 0.0: the time at which the time series begins

...

...

Outputs

  • A struct of type Demodulate

...

See also: Demodulate, and relevant plots recipe demodplot

Utilities

Multitaper.tanhtransFunction
tanhtrans(trcsq,ntf)

Inverse of the magnitude squared coherence variance stabilizing transform

...

Arguments

  • trcsq::Float64: The transformed squared coherence
  • ntf::Int64: the number of tapers use to compute the coherence

...

...

Outputs

  • The output is a Float64 indicating the reverse-transformed MSC

...

Example

If the transformed coherence is 7.3 and 6 dof, the significance is

julia> invmscsig(tanhtrans(7.3,6),6)

If the significance is 0.9 and 6 dof the transformed msc is

julia> atanhtrans(sqrt(mscsig(0.9,6)),6)

See also: multispec

Multitaper.atanhtransFunction
atanhtrans(c,ntf)

Magnitude squared coherence variance stabilizing transform

...

Arguments

  • c::Float64: the value of the coherence
  • ntf::Int64: the number of tapers use to compute the coherence

...

...

Outputs

  • The output is a Float64 indicating the transformed MSC

...

See also: multispec

Multitaper.ejnFunction
ejn(DoF)

Expected jackknife variance of an univariate multitaper spectrum estimate with DoF tapers

...

Arguments

  • DoF::Int64: the number of tapers use to compute the spectrum

...

...

Outputs

  • The output is a Float64 indicating the expected jackknife variance.

...

See also: multispec

Multitaper.unwrapphaseFunction
unwrapphase(y, typ)

Unwrap the phase

...

Arguments

  • y::Vector{Float64}: The vector of phases

  • typ::Symbol: Whether to compute the unwrapped phase in degrees or radians

...

...

Outputs

  • x::Vector{Float64}: The vector of unwrapped phases

...

Alternative phase unwrapper which takes a complex argument.

Multitaper.blockerrFunction
blockerr(lengt, nsegmentts; <keyword arguments>)

Blocker code to divide the data into segments

...

Arguments

  • lengt::Int64: input length of a time series
  • nsegments::Int64: number of segments
  • overlap::Float64 = 0.0: fraction of overlap between segments, between 0.0 and 1.0

...

...

Outputs

  • seq::Vector{Int64}: Beginning index for each segment
  • seg_len::Int64: length of the segments
  • ov::Float64: overlap that actually resulted (≈overlap, above)

...

Number of Crossings

Multitaper.Pgram_upcrossingsFunction
Pgram_upcrossings(z, N)

Number of upcrossings of the level z of a periodogram spectrum per Rayleigh resolution ...

Arguments

  • z::Float64: the level of the standardized spectrum
  • N::Int64: sampling rate in time units

...

...

Outputs

  • Float64: the number of upcrossings per Rayleigh

...

Multitaper.MT_UpcrossingsFunction
MT_Upcrossings(z, α, CR, N)

Number of upcrossings of the level z of a multitaper spectrum per Rayleigh resolution ...

Arguments

  • z::Float64: the level of the standardized spectrum
  • α::Float64: degrees of freedom of estimate
  • `CR::Float64: time bandwidth product
  • N::Int64: sampling rate in time units

...

...

Outputs

  • Float64: the number of upcrossings per Rayleigh

...

Multitaper.uctableFunction

Gives an upcrossing table for α = K, CR = NW (Float), optionally num_Ray has default 1e5, and signficiance levels that can be chosen, or set to default

Two Dimensional Spectrum

Multitaper.rectslepsFunction
rectsleps(nslep,n,m,Kp,N,M, <keyword arguments>)=)

Slepian functions concentrated in 2 dimensions on a rectangle in physical space and a circle in spectral space.

...

Arguments

  • nslep::Int64: number of output slepians
  • n::Int64: number of GL nodes in the x-direction
  • m::Int64: number of GL nodes in the y-direction
  • Kp::Float64: the radius of the circle in spectral space
  • N::Int64: number of Gauss-Legendre nodes in the first dimension
  • M::Int64: number of Gauss-Legendre nodes in the second dimension
  • verbose::Boo: select true if you would like to see the concentrations

...

...

Outputs

  • sleps::Array{Matrix{Float64},1} - an array of 2D tapers

...

Multitaper.multispec2_RectangleFunction
multispec2_Rectangle(Md, K, ntaper; <keyword arguments>)

Function for computing the 2D multitaper spectrum estimate

...

Arguments

  • Md::Matrix{T} where T<:Float64: data matrix
  • K::Float64: squared radius (bandwidth) of the region of concentration in spectrum
  • ntaper::Int64: number of tapers
  • padn::Int64 = 0: number of samples to pad to in the x-direction
  • padm::Int64 = 0: number of samples to pad to in the y-direction
  • nquad::Int64 = 72: number of quadrature nodes in the x-direction
  • mquad::Int64 = 72: number of quadrature nodes in the y-direction
  • center::Bool = true: whether to subtract the mean or not
  • sleps::Union{Matrix{T},Nothing} = nothing: if you have already computed the tapers, input them here

...

...

Outputs

  • Matrix containing the 2D multitaper spectrum (simple averaging used)

Note: If this function runs slowly for the size of your input matrix, reducing the number of quadrature nodes for computing the tapers, or alternatively precomputing the tapers may help. ...