summaryrefslogtreecommitdiff
path: root/macros/ellip.sci
blob: 1dc84356ed3844afcd77e86256a00dedebcc6aae (plain)
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
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
// Copyright (C) 2018 - IIT Bombay - FOSSEE
// This file must be used under the terms of the CeCILL.
// This source file is licensed as described in the file COPYING, which
// you should have received as part of this distribution.  The terms
// are also available at
// http://www.cecill.info/licences/Licence_CeCILL_V2-en.txt
// Original Source : https://octave.sourceforge.io/signal/
// Modifieded by: Abinash Singh Under FOSSEE Internship
// Last Modified on : 3 Feb 2024
// Organization: FOSSEE, IIT Bombay
// Email: toolbox@scilab.in
function [a, b, c, d] = ellip (n, rp, rs, w, varargin)
    //Elliptic or Cauer filter design with rp dB of passband ripple and rs dB of stopband attenuation.

    //Calling Sequence
    //[b, a] = ellip (n, rp, rs, wp)
    //[b, a] = ellip (n, rp, rs, wp, "high")
    //[b, a] = ellip (n, rp, rs, [wl, wh])
    //[b, a] = ellip (n, rp, rs, [wl, wh], "stop")
    //[z, p, g] = ellip (…)
    //[…] = ellip (…, "s")

    //Parameters
    //n: positive integer value (order of filter)
    //rp: non negative scalar value (passband ripple)
    //rs: non negative scalar value (stopband attenuation)
    //wp: positive real value,
    //    1).Normalised digital passband edge(s) for digital filter, in the range [0, 1] {dimensionless}
    //    2).Analog passband edge(s) for analog filter, in the range [0, Inf] {rad/sec}

    //Description
    //This function generates an elliptic or Cauer filter with rp dB of passband ripple and rs dB of stopband attenuation.
    //[b, a] = ellip(n, Rp, Rs, Wp) indicates low pass filter with order n, Rp decibels of ripple in the passband and a stopband Rs decibels down and cutoff of pi*Wp radians. If the fifth argument is high, then the filter is a high pass filter.
    //[b, a] = ellip(n, Rp, Rs, [Wl, Wh]) indictaes band pass filter with band pass edges pi*Wl and pi*Wh. If the fifth argument is stop, the filter is a band reject filter.
    //[z, p, g] = ellip(...) returns filter as zero-pole-gain.
    //[...] = ellip(...,’s’) returns a Laplace space filter, wp can be larger than 1.

    //Examples
    //[b, a]=ellip(2, 3, 40, [0.3,0.4])
    //Output :
    // a  =
    //
    //
    //         column 1 to 4
    //
    //    1.  - 1.7258519    2.5097172  - 1.5592802
    //
    //         column 5
    //
    //    0.8188057
    // b  =
    //
    //
    //         column 1 to 4
    //
    //    0.0202774  - 0.0164257    0.0027304  - 0.0164257
    //
    //         column 5
    //
    //    0.0202774
    // Dependencies
    // ellipap sftrans bilinear zp2tf
    funcprot(0);
    [nargout nargin] = argn();

    if (nargin > 6 | nargin < 4 | nargout > 4 | nargout < 2)
        error("ellip: invalid number of inputs");
    end

    // interpret the input parameters
    if (~ (isscalar (n) & (n == fix (n)) & (n > 0)))
        error ("ellip: filter order N must be a positive integer");
    end

    stop = %F;
    digital = %T;
    for i = 1:length(varargin)
        select (varargin(i))
        case "s"
            digital = %F;
        case "z"
            digital = %T;
        case "high"
            stop = %T;
        case "stop"
            stop = %T;
        case "low"
            stop = %F;
        case "pass"
            stop = %F;
        else
            error ("ellip: expected [high|stop] or [s|z]");
        end
    end

    [rows_w columns_w] = size(w);

    if (~ ((length (w) <= 2) & (rows_w == 1 | columns_w == 1)))
        error ("ellip: frequency must be given as WP or [WL, WH]");
    elseif ((length (w) == 2) & (w(2) <= w(1)))
        error ("ellip: W(1) must be less than W(2)");
    end

    if (digital & ~ and ((w >= 0) & (w <= 1)))
        error ("ellip: all elements of W must be in the range [0,1]");
    elseif (~ digital & ~ and (w >= 0))
        error ("ellip: all elements of W must be in the range [0,inf]");
    end

    if (~ (isscalar (rp) & or(type(rp) == [1 5 8]) & (rp >= 0)))
        error ("ellip: passband ripple RP must be a non-negative scalar");
    end

    if (~ (isscalar (rs) & or(type(rs) == [1 5 8]) & (rs >= 0)))
        error ("ellip: stopband attenuation RS must be a non-negative scalar");
    end


    // Prewarp the digital frequencies
    if (digital)
        T = 2;       // sampling frequency of 2 Hz
        w = (2 / T )* tan (%pi * w / T);
    end

    // Generate s-plane poles, zeros and gain
    [zero, pole, gain] = ellipap (n, rp, rs);
    zero = zero';
    pole = pole';

    // 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);
    end

    // convert to the correct output form
    if (nargout == 2)
        [a b] = zp2tf(zero, pole, gain);
    elseif (nargout == 3)
        a = conj(zero');
        b = conj(pole');
        c = gain;
    else
        // output ss results
        //[a, b, c, d] = zp2ss (zero, pole, gain);
        error("ellip: yet not implemented in state-space form OR invalid number of o/p arguments")
    end
endfunction

/*
%% Test input validation // all passed
error [a, b] = ellip ()
error [a, b] = ellip (1)
error [a, b] = ellip (1, 2)
error [a, b] = ellip (1, 2, 3)
error [a, b] = ellip (1, 2, 3, 4, 5, 6, 7)
error [a, b] = ellip (.5, 2, 40, .2)
error [a, b] = ellip (3, 2, 40, .2, "invalid")

test

//[b, a] = ellip (n, rp, rs, [wl, wh])
[b, a] = ellip (5, 1, 90, [.1 .2]); // passed

//[b, a] = ellip (n, rp, rs, wp)
[A, B] = ellip (6, 3, 50, .6); //passed
 
//[b, a] = ellip (n, rp, rs, [wl, wh], "stop")
[b, a] = ellip (5, 1, 90, [.1 .2],"stop"); // passed

//[b, a] = ellip (n, rp, rs, wp, "high")
[A, B] = ellip (6, 3, 50, .6,"high"); // passed

//[b, a] = ellip (n, rp, rs, [wl, wh])
[b, a] = ellip (5, 1, 90, [.1 .2],"s"); // passed

//[b, a] = ellip (n, rp, rs, wp)
[A, B] = ellip (6, 3, 50, .6,"z"); // passed
 
//[b, a] = ellip (n, rp, rs, [wl, wh], "stop")
[b, a] = ellip (5, 1, 90, [.1 .2],"pass"); // passed

//[b, a] = ellip (n, rp, rs, wp, "high")
[A, B] = ellip (6, 3, 50, .6,"low"); //  passed


//[z, p, g] = ellip (…)
[z, p, g] = ellip (6, 3, 50, .6); //  passed

[z,p,g] = ellip( 10,15,4,.9,"s") // passed

test
 [a, b, c, d] = ellip (6, 3, 50, .6); // not implemented yet 
 */