]> Creatis software - CreaPhase.git/blob - octave_packages/m/signal/periodogram.m
update packages
[CreaPhase.git] / octave_packages / m / signal / periodogram.m
1 ## Copyright (C) 1995-2012 Friedrich Leisch
2 ## Copyright (C) 2010 Alois Schloegl
3 ##
4 ## This file is part of Octave.
5 ##
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.
10 ##
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.
15 ##
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/>.
19
20 ## -*- texinfo -*-
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}.
24 ##
25 ## [Pxx,w] = periodogram (@var{x}).
26 ##
27 ## [Pxx,w] = periodogram (@var{x},win).
28 ##
29 ## [Pxx,w] = periodogram (@var{x},win,nfft).
30 ##
31 ## [Pxx,f] = periodogram (@var{x},win,nfft,Fs).
32 ##
33 ## [Pxx,f] = periodogram (@var{x},win,nfft,Fs,"range").
34 ##
35 ## @itemize
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.
39 ##
40 ## @item win: weight data with window, x.*win is used for further computation,
41 ## if window is empty, a rectangular window is used.
42 ##
43 ## @item nfft: number of frequency bins, default max(256, 2.^ceil(log2(length(x)))).
44 ##
45 ## @item Fs: sampling rate, default 1.
46 ##
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.
50 ##
51 ## @item Pxx: one-, or two-sided power spectrum.
52 ##
53 ## @item w: angular frequency [0..2*pi) (two-sided) or [0..pi] one-sided.
54 ##
55 ## @item f: frequency [0..Fs) (two-sided) or [0..Fs/2] one-sided.
56 ## @end itemize
57 ## @end deftypefn
58
59 ## Author: FL <Friedrich.Leisch@ci.tuwien.ac.at>
60 ## Description: Compute the periodogram
61
62 function [pxx, f] = periodogram (x, varargin)
63
64   ## check input arguments
65
66   if (nargin < 1 || nargin > 5)
67     print_usage ();
68   endif
69
70   nfft = []; fs = []; range = []; window = [];
71   j = 1;
72   for k = 1:length (varargin)
73     if (ischar (varargin{k}))
74       range = varargin{k};
75     else
76       switch (j)
77       case 1
78         window = varargin{k};
79       case 2
80         nfft   = varargin{k};
81       case 3
82         fs     = varargin{k};
83       case 4
84         range  = varargin{k};
85       endswitch
86       j++;
87     endif
88   endfor
89
90   [r, c] = size (x);
91   if (r == 1)
92     r = c;
93   endif
94
95   if (ischar (window))
96     range = window;
97     window = [];
98   endif;
99   if (ischar (nfft))
100     range = nfft;
101     nfft = [];
102   endif;
103   if (ischar (fs))
104     range = fs;
105     fs = [];
106   endif;
107
108   if (!  isempty (window))
109     if (all (size (x) == size (window)))
110       x .*= window;
111     elseif (size (x, 1) == size (window, 1) && size (window, 2) == 1)
112       x .*= window (:,ones (1,c));
113     endif;
114   endif
115
116   if (numel (nfft)>1)
117     error ("nfft must be scalar");
118   endif
119   if (isempty (nfft))
120     nfft = max (256, 2.^ceil (log2 (r)));
121   endif
122
123   if (strcmp (range, "onesided"))
124     range = 1;
125   elseif strcmp (range, "twosided")
126     range = 2;
127   else
128     range = 2-isreal (x);
129   endif
130
131   ## compute periodogram
132
133   if (r>nfft)
134     Pxx = 0;
135     rr = rem (length (x), nfft);
136     if (rr)
137       x = [x(:); (zeros (nfft-rr, 1))];
138     endif
139     x = sum (reshape (x, nfft, []), 2);
140   endif
141
142   if (isempty (window))
143     n = r;
144   else
145     n = sumsq (window);
146   end;
147   Pxx = (abs (fft (x, nfft))) .^ 2 / n ;
148
149   if (nargin<4)
150     Pxx /= 2*pi;
151   elseif (! isempty (fs))
152     Pxx /= fs;
153   endif
154
155   ## generate output arguments
156
157   if (range == 1)  # onesided
158     Pxx = Pxx(1:nfft/2+1) + [0; Pxx(end:-1:(nfft/2+2)); 0];
159   endif
160
161   if (nargout != 1)
162     if (range == 1)
163       f = (0:nfft/2)'/nfft;
164     elseif (range == 2)
165       f = (0:nfft-1)'/nfft;
166     endif
167     if (nargin<4)
168       f *= 2*pi; # generate w=2*pi*f
169     elseif (! isempty (fs))
170       f *= fs;
171     endif
172   endif
173
174   if (nargout == 0)
175     if (nargin<4)
176       plot (f/(2*pi), 10*log10 (Pxx));
177       xlabel ("normalized frequency [x pi rad]");
178       ylabel ("Power density [dB/rad/sample]");
179     else
180       plot (f, 10*log10 (Pxx));
181       xlabel ("frequency [Hz]");
182       ylabel ("Power density [dB/Hz]");
183     endif
184     grid on;
185     title ("Periodogram Power Spectral Density Estimate");
186   else
187     pxx = Pxx;
188   endif
189
190 endfunction