view scripts/signal/periodogram.m @ 17281:bc924baa2c4e

doc: Add new @qcode macro for code samples which are quoted. Macro handles options ("on") or properties ("position") more elegantly than @code{"text"}. * doc/interpreter/macros.texi: Add new @qcode macro. * doc/interpreter/tips.txi: Add documentation about @qcode macro. * doc/interpreter/basics.txi, doc/interpreter/container.txi, doc/interpreter/emacs.txi, doc/interpreter/errors.txi, doc/interpreter/eval.txi, doc/interpreter/expr.txi, doc/interpreter/external.txi, doc/interpreter/func.txi, doc/interpreter/grammar.txi, doc/interpreter/image.txi, doc/interpreter/install.txi, doc/interpreter/interp.txi, doc/interpreter/io.txi, doc/interpreter/matrix.txi, doc/interpreter/numbers.txi, doc/interpreter/oop.txi, doc/interpreter/package.txi, doc/interpreter/plot.txi, doc/interpreter/quad.txi, doc/interpreter/sparse.txi, doc/interpreter/strings.txi, doc/interpreter/system.txi, doc/interpreter/vectorize.txi, libinterp/corefcn/balance.cc, libinterp/corefcn/bitfcns.cc, libinterp/corefcn/cellfun.cc, libinterp/corefcn/conv2.cc, libinterp/corefcn/data.cc, libinterp/corefcn/debug.cc, libinterp/corefcn/defaults.cc, libinterp/corefcn/dirfns.cc, libinterp/corefcn/dlmread.cc, libinterp/corefcn/error.cc, libinterp/corefcn/file-io.cc, libinterp/corefcn/find.cc, libinterp/corefcn/gammainc.cc, libinterp/corefcn/graphics.cc, libinterp/corefcn/help.cc, libinterp/corefcn/hex2num.cc, libinterp/corefcn/input.cc, libinterp/corefcn/load-path.cc, libinterp/corefcn/load-save.cc, libinterp/corefcn/ls-oct-ascii.cc, libinterp/corefcn/lu.cc, libinterp/corefcn/luinc.cc, libinterp/corefcn/matrix_type.cc, libinterp/corefcn/oct-hist.cc, libinterp/corefcn/pager.cc, libinterp/corefcn/pr-output.cc, libinterp/corefcn/pt-jit.cc, libinterp/corefcn/qz.cc, libinterp/corefcn/rand.cc, libinterp/corefcn/regexp.cc, libinterp/corefcn/schur.cc, libinterp/corefcn/sighandlers.cc, libinterp/corefcn/sparse.cc, libinterp/corefcn/spparms.cc, libinterp/corefcn/str2double.cc, libinterp/corefcn/svd.cc, libinterp/corefcn/symtab.cc, libinterp/corefcn/syscalls.cc, libinterp/corefcn/toplev.cc, libinterp/corefcn/tril.cc, libinterp/corefcn/typecast.cc, libinterp/corefcn/utils.cc, libinterp/corefcn/variables.cc, libinterp/dldfcn/__init_fltk__.cc, libinterp/dldfcn/chol.cc, libinterp/dldfcn/colamd.cc, libinterp/dldfcn/fftw.cc, libinterp/dldfcn/qr.cc, libinterp/dldfcn/symbfact.cc, libinterp/octave-value/ov-base.cc, libinterp/octave-value/ov-fcn-handle.cc, libinterp/octave-value/ov-fcn-inline.cc, libinterp/octave-value/ov-java.cc, libinterp/octave-value/ov-range.cc, libinterp/octave-value/ov-struct.cc, libinterp/octave-value/ov-usr-fcn.cc, libinterp/parse-tree/oct-parse.in.yy, libinterp/parse-tree/pt-binop.cc, libinterp/parse-tree/pt-eval.cc, libinterp/parse-tree/pt-mat.cc, scripts/@ftp/ftp.m, scripts/deprecated/java_convert_matrix.m, scripts/deprecated/java_debug.m, scripts/deprecated/java_unsigned_conversion.m, scripts/deprecated/shell_cmd.m, scripts/general/dblquad.m, scripts/general/display.m, scripts/general/genvarname.m, scripts/general/idivide.m, scripts/general/interp1.m, scripts/general/interp2.m, scripts/general/interp3.m, scripts/general/interpn.m, scripts/general/isa.m, scripts/general/profexplore.m, scripts/general/profile.m, scripts/general/quadgk.m, scripts/general/randi.m, scripts/general/structfun.m, scripts/general/subsindex.m, scripts/general/triplequad.m, scripts/geometry/griddata.m, scripts/geometry/griddata3.m, scripts/geometry/griddatan.m, scripts/geometry/voronoi.m, scripts/help/help.m, scripts/help/lookfor.m, scripts/image/cmpermute.m, scripts/image/colormap.m, scripts/image/image.m, scripts/image/imagesc.m, scripts/image/imfinfo.m, scripts/image/imformats.m, scripts/image/imread.m, scripts/image/imshow.m, scripts/image/imwrite.m, scripts/image/ind2gray.m, scripts/image/lines.m, scripts/image/rgb2ind.m, scripts/image/spinmap.m, scripts/io/dlmwrite.m, scripts/io/strread.m, scripts/io/textread.m, scripts/io/textscan.m, scripts/java/javaclasspath.m, scripts/java/usejava.m, scripts/miscellaneous/bzip2.m, scripts/miscellaneous/computer.m, scripts/miscellaneous/copyfile.m, scripts/miscellaneous/debug.m, scripts/miscellaneous/dos.m, scripts/miscellaneous/edit.m, scripts/miscellaneous/gzip.m, scripts/miscellaneous/license.m, scripts/miscellaneous/mkoctfile.m, scripts/miscellaneous/movefile.m, scripts/miscellaneous/parseparams.m, scripts/miscellaneous/unix.m, scripts/optimization/fminbnd.m, scripts/optimization/fminsearch.m, scripts/optimization/fminunc.m, scripts/optimization/fsolve.m, scripts/optimization/fzero.m, scripts/optimization/glpk.m, scripts/optimization/lsqnonneg.m, scripts/optimization/optimset.m, scripts/optimization/pqpnonneg.m, scripts/pkg/pkg.m, scripts/plot/allchild.m, scripts/plot/ancestor.m, scripts/plot/area.m, scripts/plot/axis.m, scripts/plot/bar.m, scripts/plot/barh.m, scripts/plot/box.m, scripts/plot/caxis.m, scripts/plot/cla.m, scripts/plot/clabel.m, scripts/plot/clf.m, scripts/plot/close.m, scripts/plot/colorbar.m, scripts/plot/daspect.m, scripts/plot/ezmesh.m, scripts/plot/ezmeshc.m, scripts/plot/ezsurf.m, scripts/plot/ezsurfc.m, scripts/plot/findall.m, scripts/plot/findobj.m, scripts/plot/gcbo.m, scripts/plot/gcf.m, scripts/plot/gco.m, scripts/plot/grid.m, scripts/plot/guihandles.m, scripts/plot/hdl2struct.m, scripts/plot/hidden.m, scripts/plot/hold.m, scripts/plot/isonormals.m, scripts/plot/isosurface.m, scripts/plot/legend.m, scripts/plot/mesh.m, scripts/plot/meshc.m, scripts/plot/meshz.m, scripts/plot/newplot.m, scripts/plot/orient.m, scripts/plot/pareto.m, scripts/plot/patch.m, scripts/plot/pbaspect.m, scripts/plot/pcolor.m, scripts/plot/plot.m, scripts/plot/print.m, scripts/plot/private/__add_default_menu__.m, scripts/plot/quiver.m, scripts/plot/quiver3.m, scripts/plot/refreshdata.m, scripts/plot/saveas.m, scripts/plot/scatter.m, scripts/plot/scatter3.m, scripts/plot/shading.m, scripts/plot/shrinkfaces.m, scripts/plot/slice.m, scripts/plot/stem.m, scripts/plot/stem3.m, scripts/plot/struct2hdl.m, scripts/plot/subplot.m, scripts/plot/surf.m, scripts/plot/surfc.m, scripts/plot/surfl.m, scripts/plot/tetramesh.m, scripts/plot/uigetfile.m, scripts/plot/uimenu.m, scripts/plot/uiputfile.m, scripts/plot/waterfall.m, scripts/plot/whitebg.m, scripts/plot/xlim.m, scripts/plot/ylim.m, scripts/plot/zlim.m, scripts/polynomial/conv.m, scripts/polynomial/polyout.m, scripts/polynomial/splinefit.m, scripts/set/ismember.m, scripts/set/powerset.m, scripts/set/setdiff.m, scripts/set/union.m, scripts/set/unique.m, scripts/signal/detrend.m, scripts/signal/filter2.m, scripts/signal/freqz.m, scripts/signal/periodogram.m, scripts/signal/spectral_adf.m, scripts/signal/spectral_xdf.m, scripts/sparse/eigs.m, scripts/sparse/svds.m, scripts/specfun/legendre.m, scripts/special-matrix/gallery.m, scripts/statistics/base/mean.m, scripts/statistics/base/moment.m, scripts/statistics/tests/cor_test.m, scripts/statistics/tests/kolmogorov_smirnov_test.m, scripts/statistics/tests/kolmogorov_smirnov_test_2.m, scripts/statistics/tests/kruskal_wallis_test.m, scripts/statistics/tests/prop_test_2.m, scripts/statistics/tests/sign_test.m, scripts/statistics/tests/t_test.m, scripts/statistics/tests/t_test_2.m, scripts/statistics/tests/t_test_regression.m, scripts/statistics/tests/u_test.m, scripts/statistics/tests/var_test.m, scripts/statistics/tests/welch_test.m, scripts/statistics/tests/wilcoxon_test.m, scripts/statistics/tests/z_test.m, scripts/statistics/tests/z_test_2.m, scripts/strings/base2dec.m, scripts/strings/index.m, scripts/strings/isstrprop.m, scripts/strings/mat2str.m, scripts/strings/regexptranslate.m, scripts/strings/rindex.m, scripts/strings/str2num.m, scripts/strings/strcat.m, scripts/strings/strjust.m, scripts/strings/strmatch.m, scripts/strings/validatestring.m, scripts/testfun/demo.m, scripts/testfun/example.m, scripts/testfun/test.m, scripts/time/addtodate.m, scripts/time/asctime.m, scripts/time/datestr.m, scripts/time/datetick.m, scripts/time/weekday.m, scripts/ui/errordlg.m, scripts/ui/helpdlg.m, scripts/ui/inputdlg.m, scripts/ui/listdlg.m, scripts/ui/msgbox.m, scripts/ui/questdlg.m, scripts/ui/warndlg.m: Use new @qcode macro.
author Rik <rik@octave.org>
date Mon, 19 Aug 2013 20:46:38 -0700
parents c3c1ebfaa7dc
children 1c89599167a6
line wrap: on
line source

## Copyright (C) 1995-2012 Friedrich Leisch
## Copyright (C) 2010 Alois Schloegl
##
## This file is part of Octave.
##
## Octave is free software; you can redistribute it and/or modify it
## under the terms of the GNU General Public License as published by
## the Free Software Foundation; either version 3 of the License, or (at
## your option) any later version.
##
## Octave is distributed in the hope that it will be useful, but
## WITHOUT ANY WARRANTY; without even the implied warranty of
## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
## General Public License for more details.
##
## You should have received a copy of the GNU General Public License
## along with Octave; see the file COPYING.  If not, see
## <http://www.gnu.org/licenses/>.

## -*- texinfo -*-
## @deftypefn {Function File} {[Pxx, @var{w}] =} periodogram (@var{x})
## For a data matrix @var{x} from a sample of size @var{n}, return the
## periodogram.  The angular frequency is returned in @var{w}.
##
## [Pxx,w] = periodogram (@var{x}).
##
## [Pxx,w] = periodogram (@var{x},win).
##
## [Pxx,w] = periodogram (@var{x},win,nfft).
##
## [Pxx,f] = periodogram (@var{x},win,nfft,Fs).
##
## [Pxx,f] = periodogram (@var{x},win,nfft,Fs,"range").
##
## @itemize
## @item x: data; if real-valued a one-sided spectrum is estimated,
## if complex-valued or range indicates @qcode{"@nospell{twosided}"}, the full
## spectrum is estimated.
##
## @item win: weight data with window, x.*win is used for further computation,
## if window is empty, a rectangular window is used.
##
## @item nfft: number of frequency bins, default max (256, 2.^ceil (log2 (length (x)))).
##
## @item Fs: sampling rate, default 1.
##
## @item range: @qcode{"@nospell{onesided}"} computes spectrum from [0..nfft/2+1].
## @qcode{"@nospell{twosided}"} computes spectrum from [0..nfft-1].  These
## strings can appear at any position in the list input arguments after
## window.
##
## @item @nospell{Pxx}: one-, or two-sided power spectrum.
##
## @item w: angular frequency [0..2*pi) (two-sided) or [0..pi] one-sided.
##
## @item f: frequency [0..Fs) (two-sided) or [0..Fs/2] one-sided.
## @end itemize
## @end deftypefn

## Author: FL <Friedrich.Leisch@ci.tuwien.ac.at>
## Description: Compute the periodogram

function [pxx, f] = periodogram (x, varargin)

  ## check input arguments

  if (nargin < 1 || nargin > 5)
    print_usage ();
  endif

  nfft = []; fs = []; range = []; window = [];
  j = 1;
  for k = 1:length (varargin)
    if (ischar (varargin{k}))
      range = varargin{k};
    else
      switch (j)
        case 1
          window = varargin{k};
        case 2
          nfft   = varargin{k};
        case 3
          fs     = varargin{k};
        case 4
          range  = varargin{k};
      endswitch
      j++;
    endif
  endfor

  [r, c] = size (x);
  if (r == 1)
    r = c;
  endif

  if (ischar (window))
    range = window;
    window = [];
  endif;
  if (ischar (nfft))
    range = nfft;
    nfft = [];
  endif;
  if (ischar (fs))
    range = fs;
    fs = [];
  endif;

  if (!  isempty (window))
    if (all (size (x) == size (window)))
      x .*= window;
    elseif (rows (x) == rows (window) && columns (window) == 1)
      x .*= window (:,ones (1,c));
    endif;
  endif

  if (numel (nfft)>1)
    error ("nfft must be scalar");
  endif
  if (isempty (nfft))
    nfft = max (256, 2.^ceil (log2 (r)));
  endif

  if (strcmp (range, "onesided"))
    range = 1;
  elseif (strcmp (range, "twosided"))
    range = 2;
  else
    range = 2-isreal (x);
  endif

  ## compute periodogram

  if (r>nfft)
    Pxx = 0;
    rr = rem (length (x), nfft);
    if (rr)
      x = [x(:); (zeros (nfft-rr, 1))];
    endif
    x = sum (reshape (x, nfft, []), 2);
  endif

  if (isempty (window))
    n = r;
  else
    n = sumsq (window);
  end;
  Pxx = (abs (fft (x, nfft))) .^ 2 / n ;

  if (nargin<4)
    Pxx /= 2*pi;
  elseif (! isempty (fs))
    Pxx /= fs;
  endif

  ## generate output arguments

  if (range == 1)  # onesided
    Pxx = Pxx(1:nfft/2+1) + [0; Pxx(end:-1:(nfft/2+2)); 0];
  endif

  if (nargout != 1)
    if (range == 1)
      f = (0:nfft/2)'/nfft;
    elseif (range == 2)
      f = (0:nfft-1)'/nfft;
    endif
    if (nargin<4)
      f *= 2*pi; # generate w=2*pi*f
    elseif (! isempty (fs))
      f *= fs;
    endif
  endif

  if (nargout == 0)
    if (nargin<4)
      plot (f/(2*pi), 10*log10 (Pxx));
      xlabel ("normalized frequency [x pi rad]");
      ylabel ("Power density [dB/rad/sample]");
    else
      plot (f, 10*log10 (Pxx));
      xlabel ("frequency [Hz]");
      ylabel ("Power density [dB/Hz]");
    endif
    grid on;
    title ("Periodogram Power Spectral Density Estimate");
  else
    pxx = Pxx;
  endif

endfunction