Home > utilities > lanczosfilter.m

lanczosfilter

PURPOSE ^

LANCZOSFILTER Filters a time series via Lanczos method (cosine filter).

SYNOPSIS ^

function [y,coef,window,Cx,Ff] = lanczosfilter(x,varargin)

DESCRIPTION ^

LANCZOSFILTER   Filters a time series via Lanczos method (cosine filter).  
   [Y,coef,window,Cx,Ff] = LANCZOSFILTER(X,dT,Cf,M,pass) Filters the time
   series via the Lanczos filter in the frequency space (FFT), where

   INPUTS:
      X    - Time series
      dT   - Sampling interval       (default: 1)
      Cf   - Cut-off frequency       (default: half Nyquist)
      M    - Number of coefficients  (default: 10)
      pass - Low or high-pass filter (default: 'low')

   OUTPUTS:
      Y      - Filtered time series
      coef   - Coefficients of the time window (cosine)
      window - Frequency window (aprox. ones for Ff lower(greater) than Fc 
               if low(high)-pass filter and ceros otherwise)
      Cx     - Complex Fourier Transform of X for Ff frequencies
      Ff     - Fourier frequencies, from 0 to the Nyquist frequency.
  
   The program removes from the time series the frequencies greater than   
   the cut off frequency if "pass" is 'low', i.e., low-pass filter .
   Otherwise, if pass is 'high', frequencies from zero to Cf are removed,
   i.e., a high-pass filter. Units of the cut-off frequency, [Cf], must be
   [dT]^{-1}. 
  
   In consequence, when used as a low-pass the time series is smoothed   
   like a cosine filter in time space with M coefficients where greater is
   better (see the reference).    
  
   If any option is empty, defaults are used.  

   Note: NaN's elements are replaced by the mean of the time series and
         ignored. If you have a better idea, just let me know.
  
   Reference:  
   Emery, W. J. and R. E. Thomson. "Data Analysis Methods in Physical  
     Oceanography". Elsevier, 2d ed., 2004. On pages 533-539.  
  
   Example:
      dT = 30; % min 
      N = 7*24*60/dT;
      t = (0:N-1)*dT; % data for 7 days
      pnoise = 0.30;
      T1 = 12.4*60; 
      T2 = 24*60; 
      T3 = 15*24*60; 
      Tc = 10*60; % min
      xn = 5 + 3*cos(2*pi*t/T1) + 2*cos(2*pi*t/T2) + 1*cos(2*pi*t/T3);
      xn = xn + pnoise*max(xn-mean(xn))*(0.5 - rand(size(xn)));   
      [xs,c,h,Cx,f] = lanczosfilter(xn,dT,1/Tc,[],'low');  
      subplot(211), plot(t,xn,t,xs), legend('noisy','smooth'), axis tight
      subplot(212), plot(f,h,f,abs(Cx)/max(abs(Cx)),...
         [1 1]/Tc,[min(h) max(h)],'-.',...
         [1/T1 1/T2 1/T3],([1/T1 1/T2 1/T3]<=1/Tc),'o'), axis tight  
  
   See also FILTER, FFT, IFFT

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SUBFUNCTIONS ^

SOURCE CODE ^

0001 function [y,coef,window,Cx,Ff] = lanczosfilter(x,varargin)
0002 %LANCZOSFILTER   Filters a time series via Lanczos method (cosine filter).
0003 %   [Y,coef,window,Cx,Ff] = LANCZOSFILTER(X,dT,Cf,M,pass) Filters the time
0004 %   series via the Lanczos filter in the frequency space (FFT), where
0005 %
0006 %   INPUTS:
0007 %      X    - Time series
0008 %      dT   - Sampling interval       (default: 1)
0009 %      Cf   - Cut-off frequency       (default: half Nyquist)
0010 %      M    - Number of coefficients  (default: 10)
0011 %      pass - Low or high-pass filter (default: 'low')
0012 %
0013 %   OUTPUTS:
0014 %      Y      - Filtered time series
0015 %      coef   - Coefficients of the time window (cosine)
0016 %      window - Frequency window (aprox. ones for Ff lower(greater) than Fc
0017 %               if low(high)-pass filter and ceros otherwise)
0018 %      Cx     - Complex Fourier Transform of X for Ff frequencies
0019 %      Ff     - Fourier frequencies, from 0 to the Nyquist frequency.
0020 %
0021 %   The program removes from the time series the frequencies greater than
0022 %   the cut off frequency if "pass" is 'low', i.e., low-pass filter .
0023 %   Otherwise, if pass is 'high', frequencies from zero to Cf are removed,
0024 %   i.e., a high-pass filter. Units of the cut-off frequency, [Cf], must be
0025 %   [dT]^{-1}.
0026 %
0027 %   In consequence, when used as a low-pass the time series is smoothed
0028 %   like a cosine filter in time space with M coefficients where greater is
0029 %   better (see the reference).
0030 %
0031 %   If any option is empty, defaults are used.
0032 %
0033 %   Note: NaN's elements are replaced by the mean of the time series and
0034 %         ignored. If you have a better idea, just let me know.
0035 %
0036 %   Reference:
0037 %   Emery, W. J. and R. E. Thomson. "Data Analysis Methods in Physical
0038 %     Oceanography". Elsevier, 2d ed., 2004. On pages 533-539.
0039 %
0040 %   Example:
0041 %      dT = 30; % min
0042 %      N = 7*24*60/dT;
0043 %      t = (0:N-1)*dT; % data for 7 days
0044 %      pnoise = 0.30;
0045 %      T1 = 12.4*60;
0046 %      T2 = 24*60;
0047 %      T3 = 15*24*60;
0048 %      Tc = 10*60; % min
0049 %      xn = 5 + 3*cos(2*pi*t/T1) + 2*cos(2*pi*t/T2) + 1*cos(2*pi*t/T3);
0050 %      xn = xn + pnoise*max(xn-mean(xn))*(0.5 - rand(size(xn)));
0051 %      [xs,c,h,Cx,f] = lanczosfilter(xn,dT,1/Tc,[],'low');
0052 %      subplot(211), plot(t,xn,t,xs), legend('noisy','smooth'), axis tight
0053 %      subplot(212), plot(f,h,f,abs(Cx)/max(abs(Cx)),...
0054 %         [1 1]/Tc,[min(h) max(h)],'-.',...
0055 %         [1/T1 1/T2 1/T3],([1/T1 1/T2 1/T3]<=1/Tc),'o'), axis tight
0056 %
0057 %   See also FILTER, FFT, IFFT
0058 
0059 %   Written by
0060 %   Lic. on Physics Carlos Adrián Vargas Aguilera
0061 %   Physical Oceanography MS candidate
0062 %   UNIVERSIDAD DE GUADALAJARA
0063 %   Mexico, 2004
0064 %
0065 %   nubeobscura@hotmail.com
0066 
0067 % Check arguments:
0068 if nargin<1 || nargin>5
0069  error('Lanczosfilter:ArgumentNumber','Incorrect number of arguments.')
0070 elseif ~isvector(x) || ~isreal(x)
0071  error('Lanczosfilter:ArgumentType','Incorrect time series.')
0072 end
0073 if nargin<2 || isempty(varargin{1})
0074   dT = 1;
0075 elseif ~(numel(varargin{1})==1) || ~isreal(varargin{1})
0076  error('Lanczosfilter:ArgumentType','Incorrect time interval.')
0077 else
0078  dT = varargin{1};
0079 end
0080 Nf = 1/(2*dT); % Nyquist frequency
0081 if nargin<3 || isempty(varargin{2})
0082  Cf = Nf/2;
0083 elseif ~(numel(varargin{2})==1) || ~isreal(varargin{2}) || varargin{2}<=0 || varargin{2}>Nf 
0084  error('Lanczosfilter:ArgumentType','Incorrect cut-off frequency.')
0085 else
0086  Cf = varargin{2};
0087 end
0088 if nargin<4 || isempty(varargin{3})
0089  M = 100;
0090 elseif ~(numel(varargin{3})==1) || ~isreal(varargin{3}) || ~(varargin{3}==round(varargin{3}))
0091  error('Lanczosfilter:ArgumentType','Incorrect Number of coeffients.')
0092 else
0093  M = varargin{3};
0094 end
0095 if nargin<5 || isempty(varargin{4})
0096  LoH = 'l';
0097 elseif ~ischar(varargin{4}) || isempty(strfind('lh',lower(varargin{4}(1))))
0098  error('Lanczosfilter:ArgumentType','Incorrect filter pass type.')
0099 else
0100  LoH = varargin{4};
0101 end
0102 if strcmpi(LoH(1),'h')
0103  LoH = 2;
0104 else
0105  LoH = 1;
0106 end
0107 
0108 % Normalize the cut off frequency with the Nyquist frequency:
0109 Cf = Cf/Nf;
0110 
0111 % Lanczos cosine coeficients:
0112 coef = lanczos_filter_coef(Cf,M);
0113 coef = coef(:,LoH);
0114 
0115 % Filter in frequency space:
0116 [window,Ff] = spectral_window(coef,length(x)); 
0117 Ff = Ff*Nf;
0118 
0119 % Replace NaN's with the mean (ideas?):
0120 inan = isnan(x); 
0121 xmean = mean(x(~inan)); 
0122 x(inan) = xmean;
0123 
0124 % Filtering:
0125 [y,Cx] = spectral_filtering(x,window);
0126 
0127 function coef = lanczos_filter_coef(Cf,M)
0128 % Positive coeficients of Lanczos [low high]-pass.
0129 hkcs = lowpass_cosine_filter_coef(Cf,M);
0130 sigma = [1 sin(pi * (1:M) / M) ./ (pi * (1:M) / M)];
0131 hkB = hkcs.*sigma;
0132 hkA = -hkB; hkA(1) = hkA(1)+1;
0133 coef = [hkB(:) hkA(:)];
0134 
0135 function coef = lowpass_cosine_filter_coef(Cf,M)
0136 % Positive coeficients of cosine filter low-pass.
0137 coef = Cf*[1 sin(pi*(1:M)*Cf)./(pi*(1:M)*Cf)];
0138 
0139 function [window,Ff] = spectral_window(coef,N)
0140 % Window of cosine filter in frequency space.
0141 Ff = 0:2/N:1; 
0142 window = zeros(length(Ff),1);
0143 for i = 1:length(Ff)
0144  window(i) = coef(1) + 2*sum(coef(2:end).*cos((1:length(coef)-1)'*pi*Ff(i)));
0145 end
0146 
0147 function [y,Cx] = spectral_filtering(x,window)
0148 % Filtering in frequency space is multiplication, (convolution in time
0149 % space).
0150 Nx  = length(x);
0151 Cx  = fft(x(:));
0152 Cx = Cx(1:floor(Nx/2)+1);
0153 CxH = Cx.*window(:);
0154 CxH(length(CxH)+1:Nx) = conj(flipud(CxH(2:end-1)));
0155 y = real(ifft(CxH));
0156 
0157 
0158 % Carlos Adrián Vargas Aguilera. nubeobscura@hotmail.com

Generated on Wed 20-Feb-2019 16:06:01 by m2html © 2005