octave: Signal Processing

 
 31 Signal Processing
 ********************
 
 This chapter describes the signal processing and fast Fourier transform
 functions available in Octave.  Fast Fourier transforms are computed
 with the FFTW or FFTPACK libraries depending on how Octave is built.
 
  -- : fft (X)
  -- : fft (X, N)
  -- : fft (X, N, DIM)
      Compute the discrete Fourier transform of X using a Fast Fourier
      Transform (FFT) algorithm.
 
      The FFT is calculated along the first non-singleton dimension of
      the array.  Thus if X is a matrix, ‘fft (X)’ computes the FFT for
      each column of X.
 
      If called with two arguments, N is expected to be an integer
      specifying the number of elements of X to use, or an empty matrix
      to specify that its value should be ignored.  If N is larger than
      the dimension along which the FFT is calculated, then X is resized
      and padded with zeros.  Otherwise, if N is smaller than the
      dimension along which the FFT is calculated, then X is truncated.
 
      If called with three arguments, DIM is an integer specifying the
      dimension of the matrix along which the FFT is performed.
 
DONTPRINTYET       See also: Seeifft XREFifft, Seefft2 XREFfft2, *notefftn:
DONTPRINTYET       See also: Seeifft XREFifft, Seefft2 XREFfft2, Seefftn

      XREFfftn, Seefftw XREFfftw.
 
  -- : ifft (X)
  -- : ifft (X, N)
  -- : ifft (X, N, DIM)
      Compute the inverse discrete Fourier transform of X using a Fast
      Fourier Transform (FFT) algorithm.
 
      The inverse FFT is calculated along the first non-singleton
      dimension of the array.  Thus if X is a matrix, ‘fft (X)’ computes
      the inverse FFT for each column of X.
 
      If called with two arguments, N is expected to be an integer
      specifying the number of elements of X to use, or an empty matrix
      to specify that its value should be ignored.  If N is larger than
      the dimension along which the inverse FFT is calculated, then X is
      resized and padded with zeros.  Otherwise, if N is smaller than the
      dimension along which the inverse FFT is calculated, then X is
      truncated.
 
      If called with three arguments, DIM is an integer specifying the
      dimension of the matrix along which the inverse FFT is performed.
 
DONTPRINTYET       See also: Seefft XREFfft, Seeifft2 XREFifft2, *noteifftn:
DONTPRINTYET       See also: Seefft XREFfft, Seeifft2 XREFifft2, Seeifftn

      XREFifftn, Seefftw XREFfftw.
 
  -- : fft2 (A)
  -- : fft2 (A, M, N)
      Compute the two-dimensional discrete Fourier transform of A using a
      Fast Fourier Transform (FFT) algorithm.
 
      The optional arguments M and N may be used specify the number of
      rows and columns of A to use.  If either of these is larger than
      the size of A, A is resized and padded with zeros.
 
      If A is a multi-dimensional matrix, each two-dimensional sub-matrix
      of A is treated separately.
 
DONTPRINTYET       See also: Seeifft2 XREFifft2, Seefft XREFfft, *notefftn:
DONTPRINTYET       See also: Seeifft2 XREFifft2, Seefft XREFfft, Seefftn

      XREFfftn, Seefftw XREFfftw.
 
  -- : ifft2 (A)
  -- : ifft2 (A, M, N)
      Compute the inverse two-dimensional discrete Fourier transform of A
      using a Fast Fourier Transform (FFT) algorithm.
 
      The optional arguments M and N may be used specify the number of
      rows and columns of A to use.  If either of these is larger than
      the size of A, A is resized and padded with zeros.
 
      If A is a multi-dimensional matrix, each two-dimensional sub-matrix
      of A is treated separately.
 
DONTPRINTYET       See also: Seefft2 XREFfft2, Seeifft XREFifft, *noteifftn:
DONTPRINTYET       See also: Seefft2 XREFfft2, Seeifft XREFifft, Seeifftn

      XREFifftn, Seefftw XREFfftw.
 
  -- : fftn (A)
  -- : fftn (A, SIZE)
      Compute the N-dimensional discrete Fourier transform of A using a
      Fast Fourier Transform (FFT) algorithm.
 
      The optional vector argument SIZE may be used specify the
      dimensions of the array to be used.  If an element of SIZE is
      smaller than the corresponding dimension of A, then the dimension
      of A is truncated prior to performing the FFT.  Otherwise, if an
      element of SIZE is larger than the corresponding dimension then A
      is resized and padded with zeros.
 
DONTPRINTYET       See also: Seeifftn XREFifftn, Seefft XREFfft, *notefft2:
DONTPRINTYET       See also: Seeifftn XREFifftn, Seefft XREFfft, Seefft2

      XREFfft2, Seefftw XREFfftw.
 
  -- : ifftn (A)
  -- : ifftn (A, SIZE)
      Compute the inverse N-dimensional discrete Fourier transform of A
      using a Fast Fourier Transform (FFT) algorithm.
 
      The optional vector argument SIZE may be used specify the
      dimensions of the array to be used.  If an element of SIZE is
      smaller than the corresponding dimension of A, then the dimension
      of A is truncated prior to performing the inverse FFT.  Otherwise,
      if an element of SIZE is larger than the corresponding dimension
      then A is resized and padded with zeros.
 
DONTPRINTYET       See also: Seefftn XREFfftn, Seeifft XREFifft, *noteifft2:
DONTPRINTYET       See also: Seefftn XREFfftn, Seeifft XREFifft, Seeifft2

      XREFifft2, Seefftw XREFfftw.
 
    Octave uses the FFTW libraries to perform FFT computations.  When
 Octave starts up and initializes the FFTW libraries, they read a system
 wide file (on a Unix system, it is typically ‘/etc/fftw/wisdom’) that
 contains information useful to speed up FFT computations.  This
 information is called the _wisdom_.  The system-wide file allows wisdom
 to be shared between all applications using the FFTW libraries.
 
    Use the ‘fftw’ function to generate and save wisdom.  Using the
 utilities provided together with the FFTW libraries (‘fftw-wisdom’ on
 Unix systems), you can even add wisdom generated by Octave to the
 system-wide wisdom file.
 
  -- : METHOD = fftw ("planner")
  -- : fftw ("planner", METHOD)
  -- : WISDOM = fftw ("dwisdom")
  -- : fftw ("dwisdom", WISDOM)
  -- : fftw ("threads", NTHREADS)
  -- : NTHREADS = fftw ("threads")
 
      Manage FFTW wisdom data.
 
      Wisdom data can be used to significantly accelerate the calculation
      of the FFTs, but implies an initial cost in its calculation.  When
      the FFTW libraries are initialized, they read a system wide wisdom
      file (typically in ‘/etc/fftw/wisdom’), allowing wisdom to be
      shared between applications other than Octave.  Alternatively, the
      ‘fftw’ function can be used to import wisdom.  For example,
 
           WISDOM = fftw ("dwisdom")
 
      will save the existing wisdom used by Octave to the string WISDOM.
      This string can then be saved to a file and restored using the
      ‘save’ and ‘load’ commands respectively.  This existing wisdom can
      be re-imported as follows
 
           fftw ("dwisdom", WISDOM)
 
      If WISDOM is an empty string, then the wisdom used is cleared.
 
      During the calculation of Fourier transforms further wisdom is
      generated.  The fashion in which this wisdom is generated is also
      controlled by the ‘fftw’ function.  There are five different
      manners in which the wisdom can be treated:
 
      "estimate"
           Specifies that no run-time measurement of the optimal means of
           calculating a particular is performed, and a simple heuristic
           is used to pick a (probably sub-optimal) plan.  The advantage
           of this method is that there is little or no overhead in the
           generation of the plan, which is appropriate for a Fourier
           transform that will be calculated once.
 
      "measure"
           In this case a range of algorithms to perform the transform is
           considered and the best is selected based on their execution
           time.
 
      "patient"
           Similar to "measure", but a wider range of algorithms is
           considered.
 
      "exhaustive"
           Like "measure", but all possible algorithms that may be used
           to treat the transform are considered.
 
      "hybrid"
           As run-time measurement of the algorithm can be expensive,
           this is a compromise where "measure" is used for transforms up
           to the size of 8192 and beyond that the "estimate" method is
           used.
 
      The default method is "estimate".  The current method can be
      queried with
 
           METHOD = fftw ("planner")
 
      or set by using
 
           fftw ("planner", METHOD)
 
      Note that calculated wisdom will be lost when restarting Octave.
      However, the wisdom data can be reloaded if it is saved to a file
      as described above.  Saved wisdom files should not be used on
      different platforms since they will not be efficient and the point
      of calculating the wisdom is lost.
 
      The number of threads used for computing the plans and executing
      the transforms can be set with
 
           fftw ("threads", NTHREADS)
 
      Note that octave must be compiled with multi-threaded FFTW support
      for this feature.  The number of processors available to the
      current process is used per default.
 
DONTPRINTYET       See also: Seefft XREFfft, Seeifft XREFifft, *notefft2:
DONTPRINTYET DONTPRINTYET       See also: Seefft XREFfft, Seeifft XREFifft, Seefft2

      XREFfft2, Seeifft2 XREFifft2, Seefftn XREFfftn, *noteDONTPRINTYET DONTPRINTYET       See also: Seefft XREFfft, Seeifft XREFifft, Seefft2

      XREFfft2, Seeifft2 XREFifft2, Seefftn XREFfftn, See
      ifftn XREFifftn.
 
  -- : fftconv (X, Y)
  -- : fftconv (X, Y, N)
      Convolve two vectors using the FFT for computation.
 
      ‘c = fftconv (X, Y)’ returns a vector of length equal to ‘length
      (X) + length (Y) - 1’.  If X and Y are the coefficient vectors of
      two polynomials, the returned value is the coefficient vector of
      the product polynomial.
 
      The computation uses the FFT by calling the function ‘fftfilt’.  If
      the optional argument N is specified, an N-point FFT is used.
 
DONTPRINTYET       See also: Seedeconv XREFdeconv, Seeconv XREFconv, *noteDONTPRINTYET       See also: Seedeconv XREFdeconv, Seeconv XREFconv, See
      conv2 XREFconv2.
 
  -- : fftfilt (B, X)
  -- : fftfilt (B, X, N)
      Filter X with the FIR filter B using the FFT.
 
      If X is a matrix, filter each column of the matrix.
 
      Given the optional third argument, N, ‘fftfilt’ uses the
      overlap-add method to filter X with B using an N-point FFT.  The
      FFT size must be an even power of 2 and must be greater than or
      equal to the length of B.  If the specified N does not meet these
      criteria, it is automatically adjusted to the nearest value that
      does.
 
      See also: Seefilter XREFfilter, Seefilter2 XREFfilter2.
 
  -- : Y = filter (B, A, X)
  -- : [Y, SF] = filter (B, A, X, SI)
  -- : [Y, SF] = filter (B, A, X, [], DIM)
  -- : [Y, SF] = filter (B, A, X, SI, DIM)
      Apply a 1-D digital filter to the data X.
 
      ‘filter’ returns the solution to the following linear,
      time-invariant difference equation:
 
            N                   M
           SUM a(k+1) y(n-k) = SUM b(k+1) x(n-k)    for 1<=n<=length(x)
           k=0                 k=0
 
      where N=length(a)-1 and M=length(b)-1.  The result is calculated
      over the first non-singleton dimension of X or over DIM if
      supplied.
 
      An equivalent form of the equation is:
 
                     N                   M
           y(n) = - SUM c(k+1) y(n-k) + SUM d(k+1) x(n-k)  for 1<=n<=length(x)
                    k=1                 k=0
 
      where c = a/a(1) and d = b/a(1).
 
      If the fourth argument SI is provided, it is taken as the initial
      state of the system and the final state is returned as SF.  The
      state vector is a column vector whose length is equal to the length
      of the longest coefficient vector minus one.  If SI is not
      supplied, the initial state vector is set to all zeros.
 
      In terms of the Z Transform, Y is the result of passing the
      discrete-time signal X through a system characterized by the
      following rational system function:
 
                     M
                    SUM d(k+1) z^(-k)
                    k=0
           H(z) = ---------------------
                       N
                  1 + SUM c(k+1) z^(-k)
                      k=1
 
      See also: Seefilter2 XREFfilter2, Seefftfilt XREFfftfilt,
      Seefreqz XREFfreqz.
 
  -- : Y = filter2 (B, X)
  -- : Y = filter2 (B, X, SHAPE)
      Apply the 2-D FIR filter B to X.
 
      If the argument SHAPE is specified, return an array of the desired
      shape.  Possible values are:
 
      "full"
           pad X with zeros on all sides before filtering.
 
      "same"
           unpadded X (default)
 
      "valid"
           trim X after filtering so edge effects are no included.
 
      Note this is just a variation on convolution, with the parameters
      reversed and B rotated 180 degrees.
 
      See also: Seeconv2 XREFconv2.
 
  -- : [H, W] = freqz (B, A, N, "whole")
  -- : [H, W] = freqz (B)
  -- : [H, W] = freqz (B, A)
  -- : [H, W] = freqz (B, A, N)
  -- : H = freqz (B, A, W)
  -- : [H, W] = freqz (..., FS)
  -- : freqz (...)
 
      Return the complex frequency response H of the rational IIR filter
      whose numerator and denominator coefficients are B and A,
      respectively.
 
      The response is evaluated at N angular frequencies between 0 and
      2*pi.
 
      The output value W is a vector of the frequencies.
 
      If A is omitted, the denominator is assumed to be 1 (this
      corresponds to a simple FIR filter).
 
      If N is omitted, a value of 512 is assumed.  For fastest
      computation, N should factor into a small number of small primes.
 
      If the fourth argument, "whole", is omitted the response is
      evaluated at frequencies between 0 and pi.
 
      ‘freqz (B, A, W)’
 
      Evaluate the response at the specific frequencies in the vector W.
      The values for W are measured in radians.
 
      ‘[...] = freqz (..., FS)’
 
      Return frequencies in Hz instead of radians assuming a sampling
      rate FS.  If you are evaluating the response at specific
      frequencies W, those frequencies should be requested in Hz rather
      than radians.
 
      ‘freqz (...)’
 
      Plot the magnitude and phase response of H rather than returning
      them.
 
      See also: Seefreqz_plot XREFfreqz_plot.
 
  -- : freqz_plot (W, H)
  -- : freqz_plot (W, H, FREQ_NORM)
      Plot the magnitude and phase response of H.
 
      If the optional FREQ_NORM argument is true, the frequency vector W
      is in units of normalized radians.  If FREQ_NORM is false, or not
      given, then W is measured in Hertz.
 
      See also: Seefreqz XREFfreqz.
 
  -- : sinc (X)
      Compute the sinc function.
 
      Return sin (pi*x) / (pi*x).
 
  -- : B = unwrap (X)
  -- : B = unwrap (X, TOL)
  -- : B = unwrap (X, TOL, DIM)
 
      Unwrap radian phases by adding or subtracting multiples of 2*pi as
      appropriate to remove jumps greater than TOL.
 
      TOL defaults to pi.
 
      Unwrap will work along the dimension DIM.  If DIM is unspecified it
      defaults to the first non-singleton dimension.
 
  -- : [A, B] = arch_fit (Y, X, P, ITER, GAMMA, A0, B0)
      Fit an ARCH regression model to the time series Y using the scoring
      algorithm in Engle’s original ARCH paper.
 
      The model is
 
           y(t) = b(1) * x(t,1) + ... + b(k) * x(t,k) + e(t),
           h(t) = a(1) + a(2) * e(t-1)^2 + ... + a(p+1) * e(t-p)^2
 
      in which e(t) is N(0, h(t)), given a time-series vector Y up to
      time t-1 and a matrix of (ordinary) regressors X up to t.  The
      order of the regression of the residual variance is specified by P.
 
      If invoked as ‘arch_fit (Y, K, P)’ with a positive integer K, fit
      an ARCH(K, P) process, i.e., do the above with the t-th row of X
      given by
 
           [1, y(t-1), ..., y(t-k)]
 
      Optionally, one can specify the number of iterations ITER, the
      updating factor GAMMA, and initial values a0 and b0 for the scoring
      algorithm.
 
  -- : arch_rnd (A, B, T)
      Simulate an ARCH sequence of length T with AR coefficients B and CH
      coefficients A.
 
      The result y(t) follows the model
 
           y(t) = b(1) + b(2) * y(t-1) + ... + b(lb) * y(t-lb+1) + e(t),
 
      where e(t), given Y up to time t-1, is N(0, h(t)), with
 
           h(t) = a(1) + a(2) * e(t-1)^2 + ... + a(la) * e(t-la+1)^2
 
  -- : [PVAL, LM] = arch_test (Y, X, P)
      For a linear regression model
 
           y = x * b + e
 
      perform a Lagrange Multiplier (LM) test of the null hypothesis of
      no conditional heteroscedascity against the alternative of CH(P).
 
      I.e., the model is
 
           y(t) = b(1) * x(t,1) + ... + b(k) * x(t,k) + e(t),
 
      given Y up to t-1 and X up to t, e(t) is N(0, h(t)) with
 
           h(t) = v + a(1) * e(t-1)^2 + ... + a(p) * e(t-p)^2,
 
      and the null is a(1) == ... == a(p) == 0.
 
      If the second argument is a scalar integer, k, perform the same
      test in a linear autoregression model of order k, i.e., with
 
           [1, y(t-1), ..., y(t-K)]
 
      as the t-th row of X.
 
      Under the null, LM approximately has a chisquare distribution with
      P degrees of freedom and PVAL is the p-value (1 minus the CDF of
      this distribution at LM) of the test.
 
      If no output argument is given, the p-value is displayed.
 
  -- : arma_rnd (A, B, V, T, N)
      Return a simulation of the ARMA model.
 
      The ARMA model is defined by
 
           x(n) = a(1) * x(n-1) + ... + a(k) * x(n-k)
                + e(n) + b(1) * e(n-1) + ... + b(l) * e(n-l)
 
      in which K is the length of vector A, L is the length of vector B
      and E is Gaussian white noise with variance V.  The function
      returns a vector of length T.
 
      The optional parameter N gives the number of dummy X(I) used for
      initialization, i.e., a sequence of length T+N is generated and
      X(N+1:T+N) is returned.  If N is omitted, N = 100 is used.
 
  -- : autoreg_matrix (Y, K)
      Given a time series (vector) Y, return a matrix with ones in the
      first column and the first K lagged values of Y in the other
      columns.
 
      In other words, for T > K, ‘[1, Y(T-1), ..., Y(T-K)]’ is the t-th
      row of the result.
 
      The resulting matrix may be used as a regressor matrix in
      autoregressions.
 
  -- : bartlett (M)
      Return the filter coefficients of a Bartlett (triangular) window of
      length M.
 
      For a definition of the Bartlett window see, e.g., A.V. Oppenheim &
      R. W. Schafer, ‘Discrete-Time Signal Processing’.
 
  -- : blackman (M)
  -- : blackman (M, "periodic")
  -- : blackman (M, "symmetric")
      Return the filter coefficients of a Blackman window of length M.
 
      If the optional argument "periodic" is given, the periodic form of
      the window is returned.  This is equivalent to the window of length
      M+1 with the last coefficient removed.  The optional argument
      "symmetric" is equivalent to not specifying a second argument.
 
      For a definition of the Blackman window, see, e.g., A.V. Oppenheim
      & R. W. Schafer, ‘Discrete-Time Signal Processing’.
 
  -- : detrend (X, P)
      If X is a vector, ‘detrend (X, P)’ removes the best fit of a
      polynomial of order P from the data X.
 
      If X is a matrix, ‘detrend (X, P)’ does the same for each column in
      X.
 
      The second argument P is optional.  If it is not specified, a value
      of 1 is assumed.  This corresponds to removing a linear trend.
 
      The order of the polynomial can also be given as a string, in which
      case P must be either "constant" (corresponds to ‘P=0’) or "linear"
      (corresponds to ‘P=1’).
 
      See also: Seepolyfit XREFpolyfit.
 
  -- : [D, DD] = diffpara (X, A, B)
      Return the estimator D for the differencing parameter of an
      integrated time series.
 
      The frequencies from [2*pi*a/t, 2*pi*b/T] are used for the
      estimation.  If B is omitted, the interval [2*pi/T, 2*pi*a/T] is
      used.  If both B and A are omitted then a = 0.5 * sqrt (T) and b =
      1.5 * sqrt (T) is used, where T is the sample size.  If X is a
      matrix, the differencing parameter of each column is estimated.
 
      The estimators for all frequencies in the intervals described above
      is returned in DD.
 
      The value of D is simply the mean of DD.
 
      Reference: P.J. Brockwell & R.A. Davis.  ‘Time Series: Theory and
      Methods’.  Springer 1987.
 
  -- : durbinlevinson (C, OLDPHI, OLDV)
      Perform one step of the Durbin-Levinson algorithm.
 
      The vector C specifies the autocovariances ‘[gamma_0, ...,
      gamma_t]’ from lag 0 to T, OLDPHI specifies the coefficients based
      on C(T-1) and OLDV specifies the corresponding error.
 
      If OLDPHI and OLDV are omitted, all steps from 1 to T of the
      algorithm are performed.
 
  -- : fftshift (X)
  -- : fftshift (X, DIM)
      Perform a shift of the vector X, for use with the ‘fft’ and ‘ifft’
      functions, in order to move the frequency 0 to the center of the
      vector or matrix.
 
      If X is a vector of N elements corresponding to N time samples
      spaced by dt, then ‘fftshift (fft (X))’ corresponds to frequencies
 
           f = [ -(ceil((N-1)/2):-1:1), 0, (1:floor((N-1)/2)) ] * df
 
      where df = 1 / (N * dt).
 
      If X is a matrix, the same holds for rows and columns.  If X is an
      array, then the same holds along each dimension.
 
      The optional DIM argument can be used to limit the dimension along
      which the permutation occurs.
 
      See also: Seeifftshift XREFifftshift.
 
  -- : ifftshift (X)
  -- : ifftshift (X, DIM)
      Undo the action of the ‘fftshift’ function.
 
      For even length X, ‘fftshift’ is its own inverse, but odd lengths
      differ slightly.
 
      See also: Seefftshift XREFfftshift.
 
  -- : fractdiff (X, D)
      Compute the fractional differences (1-L)^d x where L denotes the
      lag-operator and d is greater than -1.
 
  -- : hamming (M)
  -- : hamming (M, "periodic")
  -- : hamming (M, "symmetric")
      Return the filter coefficients of a Hamming window of length M.
 
      If the optional argument "periodic" is given, the periodic form of
      the window is returned.  This is equivalent to the window of length
      M+1 with the last coefficient removed.  The optional argument
      "symmetric" is equivalent to not specifying a second argument.
 
      For a definition of the Hamming window see, e.g., A.V. Oppenheim &
      R. W. Schafer, ‘Discrete-Time Signal Processing’.
 
  -- : hanning (M)
  -- : hanning (M, "periodic")
  -- : hanning (M, "symmetric")
      Return the filter coefficients of a Hanning window of length M.
 
      If the optional argument "periodic" is given, the periodic form of
      the window is returned.  This is equivalent to the window of length
      M+1 with the last coefficient removed.  The optional argument
      "symmetric" is equivalent to not specifying a second argument.
 
      For a definition of the Hanning window see, e.g., A.V. Oppenheim &
      R. W. Schafer, ‘Discrete-Time Signal Processing’.
 
  -- : hurst (X)
      Estimate the Hurst parameter of sample X via the rescaled range
      statistic.
 
      If X is a matrix, the parameter is estimated for every column.
 
  -- : PP = pchip (X, Y)
  -- : YI = pchip (X, Y, XI)
      Return the Piecewise Cubic Hermite Interpolating Polynomial (pchip)
      of points X and Y.
 
      If called with two arguments, return the piecewise polynomial PP
      that may be used with ‘ppval’ to evaluate the polynomial at
      specific points.
 
      When called with a third input argument, ‘pchip’ evaluates the
      pchip polynomial at the points XI.  The third calling form is
      equivalent to ‘ppval (pchip (X, Y), XI)’.
 
      The variable X must be a strictly monotonic vector (either
      increasing or decreasing) of length N.
 
      Y can be either a vector or array.  If Y is a vector then it must
      be the same length N as X.  If Y is an array then the size of Y
      must have the form ‘[S1, S2, ..., SK, N]’ The array is reshaped
      internally to a matrix where the leading dimension is given by ‘S1
      * S2 * ... * SK’ and each row of this matrix is then treated
      separately.  Note that this is exactly opposite to ‘interp1’ but is
      done for MATLAB compatibility.
 
DONTPRINTYET       See also: Seespline XREFspline, Seeppval XREFppval, *noteDONTPRINTYET       See also: Seespline XREFspline, Seeppval XREFppval, See
      mkpp XREFmkpp, Seeunmkpp XREFunmkpp.
 
  -- : [PXX, W] = periodogram (X)
  -- : [PXX, W] = periodogram (X, WIN)
  -- : [PXX, W] = periodogram (X, WIN, NFFT)
  -- : [PXX, F] = periodogram (X, WIN, NFFT, FS)
  -- : [PXX, F] = periodogram (..., "RANGE")
  -- : periodogram (...)
      Return the periodogram (Power Spectral Density) of X.
 
      The possible inputs are:
 
      X
 
           data vector.  If X is real-valued a one-sided spectrum is
           estimated.  If X is complex-valued, or "RANGE" specifies
           "twosided", the full spectrum is estimated.
 
      WIN
           window weight data.  If window is empty or unspecified a
           default rectangular window is used.  Otherwise, the window is
           applied to the signal (‘X .* WIN’) before computing the
           periodogram.  The window data must be a vector of the same
           length as X.
 
      NFFT
           number of frequency bins.  The default is 256 or the next
           higher power of 2 greater than the length of X (‘max (256,
           2.^nextpow2 (length (x)))’).  If NFFT is greater than the
           length of the input then X will be zero-padded to the length
           of NFFT.
 
      FS
           sampling rate.  The default is 1.
 
      RANGE
           range of spectrum.  "onesided" computes spectrum from
           [0:nfft/2+1].  "twosided" computes spectrum from [0:nfft-1].
 
      The optional second output W are the normalized angular
      frequencies.  For a one-sided calculation W is in the range [0, pi]
      if NFFT is even and [0, pi) if NFFT is odd.  Similarly, for a
      two-sided calculation W is in the range [0, 2*pi] or [0, 2*pi)
      depending on NFFT.
 
      If a sampling frequency is specified, FS, then the output
      frequencies F will be in the range [0, FS/2] or [0, FS/2) for
      one-sided calculations.  For two-sided calculations the range will
      be [0, FS).
 
      When called with no outputs the periodogram is immediately plotted
      in the current figure window.
 
      See also: Seefft XREFfft.
 
  -- : sinetone (FREQ, RATE, SEC, AMPL)
      Return a sinetone of frequency FREQ with a length of SEC seconds at
      sampling rate RATE and with amplitude AMPL.
 
      The arguments FREQ and AMPL may be vectors of common size.
 
      The defaults are RATE = 8000, SEC = 1, and AMPL = 64.
 
      See also: Seesinewave XREFsinewave.
 
  -- : sinewave (M, N, D)
      Return an M-element vector with I-th element given by ‘sin (2 * pi
      * (I+D-1) / N)’.
 
      The default value for D is 0 and the default value for N is M.
 
      See also: Seesinetone XREFsinetone.
 
  -- : spectral_adf (C)
  -- : spectral_adf (C, WIN)
  -- : spectral_adf (C, WIN, B)
      Return the spectral density estimator given a vector of
      autocovariances C, window name WIN, and bandwidth, B.
 
      The window name, e.g., "triangle" or "rectangle" is used to search
      for a function called ‘WIN_lw’.
 
      If WIN is omitted, the triangle window is used.
 
      If B is omitted, ‘1 / sqrt (length (X))’ is used.
 
      See also: Seespectral_xdf XREFspectral_xdf.
 
  -- : spectral_xdf (X)
  -- : spectral_xdf (X, WIN)
  -- : spectral_xdf (X, WIN, B)
      Return the spectral density estimator given a data vector X, window
      name WIN, and bandwidth, B.
 
      The window name, e.g., "triangle" or "rectangle" is used to search
      for a function called ‘WIN_sw’.
 
      If WIN is omitted, the triangle window is used.
 
      If B is omitted, ‘1 / sqrt (length (X))’ is used.
 
      See also: Seespectral_adf XREFspectral_adf.
 
  -- : spencer (X)
      Return Spencer’s 15 point moving average of each column of X.
 
  -- : Y = stft (X)
  -- : Y = stft (X, WIN_SIZE)
  -- : Y = stft (X, WIN_SIZE, INC)
  -- : Y = stft (X, WIN_SIZE, INC, NUM_COEF)
  -- : Y = stft (X, WIN_SIZE, INC, NUM_COEF, WIN_TYPE)
  -- : [Y, C] = stft (...)
      Compute the short-time Fourier transform of the vector X with
      NUM_COEF coefficients by applying a window of WIN_SIZE data points
      and an increment of INC points.
 
      Before computing the Fourier transform, one of the following
      windows is applied:
 
      "hanning"
           win_type = 1
 
      "hamming"
           win_type = 2
 
      "rectangle"
           win_type = 3
 
      The window names can be passed as strings or by the WIN_TYPE
      number.
 
      The following defaults are used for unspecified arguments: WIN_SIZE
      = 80, INC = 24, NUM_COEF = 64, and WIN_TYPE = 1.
 
      ‘Y = stft (X, ...)’ returns the absolute values of the Fourier
      coefficients according to the NUM_COEF positive frequencies.
 
      ‘[Y, C] = stft (x, ...)’ returns the entire STFT-matrix Y and a
      3-element vector C containing the window size, increment, and
      window type, which is needed by the ‘synthesis’ function.
 
      See also: Seesynthesis XREFsynthesis.
 
  -- : X = synthesis (Y, C)
      Compute a signal from its short-time Fourier transform Y and a
      3-element vector C specifying window size, increment, and window
      type.
 
      The values Y and C can be derived by
 
           [Y, C] = stft (X , ...)
 
      See also: Seestft XREFstft.
 
  -- : [A, V] = yulewalker (C)
      Fit an AR (p)-model with Yule-Walker estimates given a vector C of
      autocovariances ‘[gamma_0, ..., gamma_p]’.
 
      Returns the AR coefficients, A, and the variance of white noise, V.