This file is indexed.

/usr/share/octave/packages/signal-1.3.0/ellip.m is in octave-signal 1.3.0-1.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
## Copyright (C) 2001 Paulo Neis <p_neis@yahoo.com.br>
## Copyright (C) 2003 Doug Stewart <dastew@sympatico.ca>
##
## This program 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.
##
## This program 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
## this program; if not, see <http://www.gnu.org/licenses/>.

## -*- texinfo -*-
## @deftypefn  {Function File} {[@var{b}, @var{a}] =} ellip (@var{n}, @var{Rp}, @var{Rs}, @var{Wp})
## @deftypefnx {Function File} {[@var{b}, @var{a}] =} ellip (@var{n}, @var{Rp}, @var{Rs}, @var{Wp}, "high")
## @deftypefnx {Function File} {[@var{b}, @var{a}] =} ellip (@var{n}, @var{Rp}, @var{Rs}, @var{[Wl}, @var{Wh}])
## @deftypefnx {Function File} {[@var{b}, @var{a}] =} ellip (@var{n}, @var{Rp}, @var{Rs}, @var{[Wl}, @var{Wh}], "stop")
## @deftypefnx {Function File} {[@var{z}, @var{p}, @var{g}] =} ellip (@dots{})
## @deftypefnx {Function File} {[@var{a}, @var{b}, @var{c}, @var{d}] =} ellip (@dots{})
## @deftypefnx {Function File} {[@dots{}] =} ellip (@dots{}, "s")
##
## Generate an Elliptic or Cauer filter (discrete and contnuious).
##
## [b,a] = ellip(n, Rp, Rs, Wp)
##  low pass filter with order n, cutoff pi*Wp radians, Rp decibels
##  of ripple in the passband and a stopband Rs decibels down.
##
## [b,a] = ellip(n, Rp, Rs, Wp, 'high')
##  high pass filter with cutoff pi*Wp...
##
## [b,a] = ellip(n, Rp, Rs, [Wl, Wh])
##  band pass filter with band pass edges pi*Wl and pi*Wh ...
##
## [b,a] = ellip(n, Rp, Rs, [Wl, Wh], 'stop')
##  band reject filter with edges pi*Wl and pi*Wh, ...
##
## [z,p,g] = ellip(...)
##  return filter as zero-pole-gain.
##
## [...] = ellip(...,'s')
##     return a Laplace space filter, W can be larger than 1.
##
## [a,b,c,d] = ellip(...)
##  return  state-space matrices
##
## References:
##
## - Oppenheim, Alan V., Discrete Time Signal Processing, Hardcover, 1999.
## - Parente Ribeiro, E., Notas de aula da disciplina TE498 -  Processamento
##   Digital de Sinais, UFPR, 2001/2002.
## - Kienzle, Paul, functions from Octave-Forge, 1999 (http://octave.sf.net).
## @end deftypefn

function [a,b,c,d] = ellip(n, Rp, Rs, W, varargin)

  if (nargin>6 || nargin<4) || (nargout>4 || nargout<2)
    print_usage;
  endif

  ## interpret the input parameters
  if (!(length(n)==1 && n == round(n) && n > 0))
    error ("ellip: filter order n must be a positive integer");
  endif


  stop = 0;
  digital = 1;
  for i=1:length(varargin)
    switch varargin{i}
      case 's', digital = 0;
      case 'z', digital = 1;
      case { 'high', 'stop' }, stop = 1;
      case { 'low',  'pass' }, stop = 0;
      otherwise,  error ("ellip: expected [high|stop] or [s|z]");
    endswitch
  endfor

  [r, c]=size(W);
  if (!(length(W)<=2 && (r==1 || c==1)))
    error ("ellip: frequency must be given as w0 or [w0, w1]");
  elseif (!(length(W)==1 || length(W) == 2))
    error ("ellip: only one filter band allowed");
  elseif (length(W)==2 && !(W(1) < W(2)))
    error ("ellip: first band edge must be smaller than second");
  endif

  if ( digital && !all(W >= 0 & W <= 1))
    error ("ellip: critical frequencies must be in (0 1)");
  elseif ( !digital && !all(W >= 0 ))
    error ("ellip: critical frequencies must be in (0 inf)");
  endif

  if (Rp < 0)
    error("ellip: passband ripple must be positive decibels");
  endif

  if (Rs < 0)
    error("ellip: stopband ripple must be positive decibels");
  endif


  ## Prewarp the digital frequencies
  if digital
    T = 2;       # sampling frequency of 2 Hz
    W = tan(pi*W/T);
  endif

  ## Generate s-plane poles, zeros and gain
  [zero, pole, gain] = ncauer(Rp, Rs, n);

  ## splane frequency transform
  [zero, pole, gain] = sftrans(zero, pole, gain, W, stop);

  ## Use bilinear transform to convert poles to the z plane
  if digital
    [zero, pole, gain] = bilinear(zero, pole, gain, T);
  endif


  ## convert to the correct output form
  if nargout==2,
    a = real(gain*poly(zero));
    b = real(poly(pole));
  elseif nargout==3,
    a = zero;
    b = pole;
    c = gain;
  else
    ## output ss results
    [a, b, c, d] = zp2ss (zero, pole, gain);
  endif

endfunction

%!demo
%! [n, Ws] = ellipord ([.1 .2], .4, 1, 90);
%! [b, a] = ellip (5, 1, 90, [.1 .2]);
%! [h, w] = freqz (b, a);
%!
%! plot (w./pi, 20*log10 (abs (h)), ";;")
%! xlabel ("Frequency");
%! ylabel ("abs(H[w])[dB]");
%! axis ([0, 1, -100, 0]);
%!
%! hold ("on");
%! x=ones (1, length (h));
%! plot (w./pi, x.*-1, ";-1 dB;")
%! plot (w./pi, x.*-90, ";-90 dB;")
%! hold ("off");