1 ## Copyright (C) 1995-2012 Friedrich Leisch
2 ## Copyright (C) 2010 Alois Schloegl
4 ## This file is part of Octave.
6 ## Octave is free software; you can redistribute it and/or modify it
7 ## under the terms of the GNU General Public License as published by
8 ## the Free Software Foundation; either version 3 of the License, or (at
9 ## your option) any later version.
11 ## Octave is distributed in the hope that it will be useful, but
12 ## WITHOUT ANY WARRANTY; without even the implied warranty of
13 ## MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 ## General Public License for more details.
16 ## You should have received a copy of the GNU General Public License
17 ## along with Octave; see the file COPYING. If not, see
18 ## <http://www.gnu.org/licenses/>.
21 ## @deftypefn {Function File} {[Pxx, @var{w}] =} periodogram (@var{x})
22 ## For a data matrix @var{x} from a sample of size @var{n}, return the
23 ## periodogram. The angular frequency is returned in @var{w}.
25 ## [Pxx,w] = periodogram (@var{x}).
27 ## [Pxx,w] = periodogram (@var{x},win).
29 ## [Pxx,w] = periodogram (@var{x},win,nfft).
31 ## [Pxx,f] = periodogram (@var{x},win,nfft,Fs).
33 ## [Pxx,f] = periodogram (@var{x},win,nfft,Fs,"range").
36 ## @item x: data; if real-valued a one-sided spectrum is estimated,
37 ## if complex-valued or range indicates "@nospell{twosided}", the full
38 ## spectrum is estimated.
40 ## @item win: weight data with window, x.*win is used for further computation,
41 ## if window is empty, a rectangular window is used.
43 ## @item nfft: number of frequency bins, default max(256, 2.^ceil(log2(length(x)))).
45 ## @item Fs: sampling rate, default 1.
47 ## @item range: "@nospell{onesided}" computes spectrum from [0..nfft/2+1].
48 ## "@nospell{twosided}" computes spectrum from [0..nfft-1]. These strings
49 ## can appear at any position in the list input arguments after window.
51 ## @item Pxx: one-, or two-sided power spectrum.
53 ## @item w: angular frequency [0..2*pi) (two-sided) or [0..pi] one-sided.
55 ## @item f: frequency [0..Fs) (two-sided) or [0..Fs/2] one-sided.
59 ## Author: FL <Friedrich.Leisch@ci.tuwien.ac.at>
60 ## Description: Compute the periodogram
62 function [pxx, f] = periodogram (x, varargin)
64 ## check input arguments
66 if (nargin < 1 || nargin > 5)
70 nfft = []; fs = []; range = []; window = [];
72 for k = 1:length (varargin)
73 if (ischar (varargin{k}))
108 if (! isempty (window))
109 if (all (size (x) == size (window)))
111 elseif (size (x, 1) == size (window, 1) && size (window, 2) == 1)
112 x .*= window (:,ones (1,c));
117 error ("nfft must be scalar");
120 nfft = max (256, 2.^ceil (log2 (r)));
123 if (strcmp (range, "onesided"))
125 elseif strcmp (range, "twosided")
128 range = 2-isreal (x);
131 ## compute periodogram
135 rr = rem (length (x), nfft);
137 x = [x(:); (zeros (nfft-rr, 1))];
139 x = sum (reshape (x, nfft, []), 2);
142 if (isempty (window))
147 Pxx = (abs (fft (x, nfft))) .^ 2 / n ;
151 elseif (! isempty (fs))
155 ## generate output arguments
157 if (range == 1) # onesided
158 Pxx = Pxx(1:nfft/2+1) + [0; Pxx(end:-1:(nfft/2+2)); 0];
163 f = (0:nfft/2)'/nfft;
165 f = (0:nfft-1)'/nfft;
168 f *= 2*pi; # generate w=2*pi*f
169 elseif (! isempty (fs))
176 plot (f/(2*pi), 10*log10 (Pxx));
177 xlabel ("normalized frequency [x pi rad]");
178 ylabel ("Power density [dB/rad/sample]");
180 plot (f, 10*log10 (Pxx));
181 xlabel ("frequency [Hz]");
182 ylabel ("Power density [dB/Hz]");
185 title ("Periodogram Power Spectral Density Estimate");