File Coverage

blib/lib/PDL/DSP/Windows.pm
Criterion Covered Total %
statement 519 560 92.6
branch 229 266 86.0
condition 32 41 78.0
subroutine 108 112 96.4
pod 58 97 59.7
total 946 1076 87.9


line stmt bran cond sub pod time code
1             package PDL::DSP::Windows;
2              
3             our $VERSION = '0.102';
4              
5 7     7   1269480 use strict;
  7         13  
  7         301  
6 7     7   35 use warnings;
  7         39  
  7         341  
7              
8 7     7   2972 use PDL::Bad ();
  7         335668  
  7         234  
9 7     7   3293 use PDL::Basic ();
  7         47863  
  7         215  
10 7     7   61 use PDL::Core ();
  7         18  
  7         160  
11 7     7   37 use PDL::Math ();
  7         13  
  7         153  
12 7     7   38 use PDL::MatrixOps ();
  7         22  
  7         177  
13 7     7   4675 use PDL::Ops ();
  7         58136  
  7         240  
14 7     7   52 use PDL::Options ();
  7         13  
  7         127  
15 7     7   38 use PDL::Primitive ();
  7         13  
  7         125  
16 7     7   4253 use PDL::Ufunc ();
  7         26931  
  7         471  
17              
18             # These constants are deleted at the end of this package
19             use constant {
20 7         36 I => PDL::Core::pdl('i'),
21 7     7   58 };
  7         14  
22              
23             # These constants are left in our namespace for historical reasons
24 7     7   5548 use constant PI => 4 * atan2(1, 1);
  7         16  
  7         504  
25 7     7   93 use constant TPI => 2 * PI;
  7         16  
  7         430  
26              
27 7     7   128 use Exporter 'import';
  7         29  
  7         89262  
28              
29             my %symmetric_windows = (
30             bartlett => 1,
31             bartlett_hann => 1,
32             blackman => 1,
33             blackman_bnh => 1,
34             blackman_ex => 1,
35             blackman_gen => 1,
36             blackman_gen3 => 1,
37             blackman_gen4 => 1,
38             blackman_gen5 => 1,
39             blackman_harris => 1,
40             blackman_harris4 => 1,
41             blackman_nuttall => 1,
42             bohman => 1,
43             cauchy => 1,
44             chebyshev => 1,
45             cos_alpha => 1,
46             cosine => 1,
47             dpss => 1,
48             exponential => 1,
49             flattop => 1,
50             gaussian => 1,
51             hamming => 1,
52             hamming_ex => 1,
53             hamming_gen => 1,
54             hann => 1,
55             hann_matlab => 1,
56             hann_poisson => 1,
57             kaiser => 1,
58             lanczos => 1,
59             nuttall => 1,
60             nuttall1 => 1,
61             parzen => 1,
62             parzen_octave => 1,
63             poisson => 1,
64             rectangular => 1,
65             triangular => 1,
66             tukey => 1,
67             welch => 1,
68             );
69              
70             # These are not exported
71             my %periodic_windows = (
72             bartlett => 1,
73             bartlett_hann => 1,
74             blackman => 1,
75             blackman_bnh => 1,
76             blackman_ex => 1,
77             blackman_gen => 1,
78             blackman_gen3 => 1,
79             blackman_gen4 => 1,
80             blackman_gen5 => 1,
81             blackman_harris => 1,
82             blackman_harris4 => 1,
83             blackman_nuttall => 1,
84             bohman => 1,
85             cauchy => 1,
86             cos_alpha => 1,
87             cosine => 1,
88             dpss => 1,
89             exponential => 1,
90             flattop => 1,
91             gaussian => 1,
92             hamming => 1,
93             hamming_ex => 1,
94             hamming_gen => 1,
95             hann => 1,
96             hann_poisson => 1,
97             kaiser => 1,
98             lanczos => 1,
99             nuttall => 1,
100             nuttall1 => 1,
101             parzen => 1,
102             poisson => 1,
103             rectangular => 1,
104             triangular => 1,
105             tukey => 1,
106             welch => 1,
107             );
108              
109             my %window_aliases = (
110             bartlett_hann => [ 'Modified Bartlett-Hann' ],
111             bartlett => [ 'fejer' ],
112             blackman_harris4 => [ 'Blackman-Harris' ],
113             blackman_harris => [ 'Minimum three term (sample) Blackman-Harris' ],
114             cauchy => [ 'Abel', 'Poisson' ],
115             chebyshev => [ 'Dolph-Chebyshev' ],
116             cos_alpha => [ 'Power-of-cosine' ],
117             cosine => [ 'sine' ],
118             dpss => [ 'sleppian' ],
119             gaussian => [ 'Weierstrass' ],
120             hann => [ 'hanning' ],
121             kaiser => [ 'Kaiser-Bessel' ],
122             lanczos => [ 'sinc' ],
123             parzen => [ 'Jackson', 'Valle-Poussin' ],
124             rectangular => [ 'dirichlet', 'boxcar' ],
125             tukey => [ 'tapered cosine' ],
126             welch => [ 'Riez', 'Bochner', 'Parzen', 'parabolic' ],
127             );
128              
129             my %window_parameters = (
130             blackman_gen => [ '$alpha' ],
131             blackman_gen3 => [ '$a0', '$a1', '$a2' ],
132             blackman_gen4 => [ '$a0', '$a1', '$a2', '$a3' ],
133             blackman_gen5 => [ '$a0', '$a1', '$a2', '$a3', '$a4' ],
134             cauchy => [ '$alpha' ],
135             chebyshev => [ '$at' ],
136             cos_alpha => [ '$alpha' ],
137             dpss => [ '$beta' ],
138             gaussian => [ '$beta' ],
139             hamming_gen => [ '$a' ],
140             hann_poisson => [ '$alpha' ],
141             kaiser => [ '$beta' ],
142             poisson => [ '$alpha' ],
143             tukey => [ '$alpha' ],
144             );
145              
146             my %window_names = (
147             bartlett_hann => 'Bartlett-Hann',
148             blackman_bnh => '*An improved version of the 3-term Blackman-Harris window given by Nuttall (Ref 2, p. 89).',
149             blackman_ex => q!'exact' Blackman!,
150             blackman_gen => '*A single parameter family of the 3-term Blackman window. ',
151             blackman_gen3 => '*The general form of the Blackman family. ',
152             blackman_gen4 => '*The general 4-term Blackman-Harris window. ',
153             blackman_gen5 => '*The general 5-term Blackman-Harris window. ',
154             blackman_harris4 => 'minimum (sidelobe) four term Blackman-Harris',
155             blackman_harris => 'Blackman-Harris',
156             blackman_nuttall => 'Blackman-Nuttall',
157             blackman => q!'classic' Blackman!,
158             dpss => 'Digital Prolate Spheroidal Sequence (DPSS)',
159             flattop => 'flat top',
160             hamming_ex => q!'exact' Hamming!,
161             hamming_gen => 'general Hamming',
162             hann_matlab => '*Equivalent to the Hann window of N+2 points, with the endpoints (which are both zero) removed.',
163             hann_poisson => 'Hann-Poisson',
164             nuttall1 => '*A window referred to as the Nuttall window.',
165             parzen_octave => 'Parzen',
166             );
167              
168             my %window_print_names = (
169             blackman_bnh => 'Blackman-Harris (bnh)',
170             blackman_gen => 'General classic Blackman',
171             hann_matlab => 'Hann (matlab)',
172             nuttall1 => 'Nuttall (v1)',
173             );
174              
175             our @EXPORT_OK = (
176             keys %symmetric_windows,
177             qw( window list_windows chebpoly cos_mult_to_pow cos_pow_to_mult ),
178             );
179              
180             $PDL::onlinedoc->scan(__FILE__) if $PDL::onlinedoc;
181              
182             =encoding utf8
183              
184             =head1 NAME
185              
186             PDL::DSP::Windows - Window functions for signal processing
187              
188             =head1 SYNOPSIS
189              
190             use PDL;
191             use PDL::DSP::Windows 'window';
192              
193             # Get an ndarray with a window's samples with the helper
194             my $samples = window( 10, tukey => { params => .5 });
195              
196             # Or construct a window object with the same parameters
197             my $window = PDL::DSP::Windows->new( 10, tukey => { params => .5 });
198              
199             # These two are equivalent
200             $samples = $window->samples;
201              
202             # The window object gives access to additional methods
203             print $window->coherent_gain, "\n";
204              
205             $window->plot; # Requires PDL::Graphics::Simple
206              
207             =head1 DESCRIPTION
208              
209             This module provides symmetric and periodic (DFT-symmetric) window functions
210             for use in filtering and spectral analysis. It provides a high-level access
211             subroutine L. This functional interface is sufficient for getting
212             the window samples. For analysis and plotting, etc. an object oriented
213             interface is provided. The functional subroutines must be either explicitly
214             exported, or fully qualified. In this document, the word I refers
215             only to the mathematical window functions, while the word I is
216             used to describe code.
217              
218             Window functions are also known as apodization functions or tapering functions.
219             In this module, each of these functions maps a sequence of C<$N> integers to
220             values called a B. (To confuse matters, the word I also has
221             other meanings when describing window functions.) The functions are often
222             named for authors of journal articles. Be aware that across the literature
223             and software, some functions referred to by several different names, and some
224             names refer to several different functions. As a result, the choice of window
225             names is somewhat arbitrary.
226              
227             The L window function requires L. The
228             L window function requires L. But the
229             remaining window functions may be used if these modules are not installed.
230              
231             The most common and easiest usage of this module is indirect, via some
232             higher-level filtering interface, such as L. The next
233             easiest usage is to return a pdl of real-space samples with the subroutine
234             L. Finally, for analyzing window functions, object methods, such as
235             L, L, L are provided.
236              
237             In the following, first the functional interface (non-object oriented) is
238             described in L. Next, the object methods are
239             described in L. Next the low-level subroutines returning samples
240             for each named window are described in L. Finally, some
241             support routines that may be of interest are described in
242             L.
243              
244             =head1 FUNCTIONAL INTERFACE
245              
246             =head2 window
247              
248             $win = window({ OPTIONS });
249             $win = window( $N, { OPTIONS });
250             $win = window( $N, $name, { OPTIONS });
251             $win = window( $N, $name, $params, { OPTIONS });
252             $win = window( $N, $name, $params, $periodic );
253              
254             Returns an C<$N> point window of type C<$name>. The arguments may be passed
255             positionally in the order C<$N, $name, $params, $periodic>, or they may be
256             passed by name in the hash C.
257              
258             =head3 EXAMPLES
259              
260             # Each of the following return a 100 point symmetric hamming window.
261              
262             $win = window(100);
263             $win = window( 100, 'hamming' );
264             $win = window( 100, { name => 'hamming' });
265             $win = window({ N => 100, name => 'hamming' });
266              
267             # Each of the following returns a 100 point symmetric hann window.
268              
269             $win = window( 100, 'hann' );
270             $win = window( 100, { name => 'hann' });
271              
272             # Returns a 100 point periodic hann window.
273              
274             $win = window( 100, 'hann', { periodic => 1 });
275              
276             # Returns a 100 point symmetric Kaiser window with alpha = 2.
277              
278             $win = window( 100, 'kaiser', { params => 2 });
279              
280             =head3 OPTIONS
281              
282             The options follow default PDL::Options rules. They may be abbreviated, and
283             are case-insensitive.
284              
285             =over
286              
287             =item B
288              
289             (string) name of window function. Default: C. This selects one of
290             the window functions listed below. Note that the suffix '_per', for periodic,
291             may be ommitted. It is specified with the option C<< periodic => 1 >>
292              
293             =item B
294              
295             ref to array of parameter or parameters for the window-function subroutine.
296             Only some window-function subroutines take parameters. If the subroutine takes
297             a single parameter, it may be given either as a number, or a list of one
298             number. For example C<3> or C<[3]>.
299              
300             =item B
301              
302             number of points in window function (the same as the order of the filter).
303             As of 0.101, throws exception if the value for C is undefined or zero.
304              
305             =item B
306              
307             If value is true, return a periodic rather than a symmetric window function.
308             Defaults to false, meaning "symmetric".
309              
310             =back
311              
312             =cut
313              
314 167     167 1 592623 sub window { PDL::DSP::Windows->new(@_)->samples }
315              
316             =head2 list_windows
317              
318             print join ", ", list_windows(), "\n"
319             print join ", ", list_windows(STR), "\n"
320              
321             C returns the names all of the available windows.
322             C returns only the names of windows matching the
323             regular expression C.
324              
325             =cut
326              
327             sub list_windows {
328 7     7 1 288039 my ($expr) = @_;
329 7 100       37 return sort keys %symmetric_windows if !$expr;
330 6         7 my @match;
331 6         127 for my $name ( sort keys %symmetric_windows ) {
332 228 100       522 if ( $name =~ /$expr/ ) {
333 59         70 push @match, $name;
334 59         56 next;
335             }
336             push @match,
337             map "$name (alias $_)",
338 169   100     227 grep /$expr/i, @{ $window_aliases{$name} // [] };
  169         502  
339             }
340 6         44 @match;
341             }
342              
343             =head1 METHODS
344              
345             =head2 new
346              
347             =for usage
348              
349             my $win = PDL::DSP::Windows->new(ARGS);
350              
351             =for ref
352              
353             Create an instance of a window object. If C are given, the instance
354             is initialized. C are interpreted in exactly the same way as arguments
355             for the L subroutine.
356              
357             =for example
358              
359             For example:
360              
361             my $win1 = PDL::DSP::Windows->new( 8, 'hann' );
362             my $win2 = PDL::DSP::Windows->new({ N => 8, name => 'hann' });
363              
364             =cut
365              
366             sub new {
367 203     203 1 1048644 my $proto = shift;
368 203   66     1089 my $self = bless {}, ref $proto || $proto;
369 203 100       851 $self->init(@_) if @_;
370 202         667 return $self;
371             }
372              
373             =head2 init
374              
375             =for usage
376              
377             $win->init(ARGS);
378              
379             =for ref
380              
381             Initialize (or reinitialize) a window object. C are interpreted in
382             exactly the same way as arguments for the L subroutine.
383             As of 0.101, throws exception if the value for C is undefined or zero.
384              
385             =for example
386              
387             For example:
388              
389             my $win = PDL::DSP::Windows->new( 8, 'hann' );
390             $win->init( 10, 'hamming' );
391              
392             =cut
393              
394             sub init {
395 252     252 1 77236 my $self = shift;
396              
397 252         438 my ( $N, $name, $params, $periodic );
398              
399 252 100       620 $N = shift unless ref $_[0];
400 252 100       547 $name = shift unless ref $_[0];
401 252 100       744 $params = shift unless ref $_[0] eq 'HASH';
402 252 100       556 $periodic = shift unless ref $_[0];
403              
404 252   100     1906 my $opts = PDL::Options->new({
405             name => 'hamming',
406             periodic => 0, # symmetric or periodic
407             N => undef, # order
408             params => undef,
409             })->options( shift // {} );
410              
411 252   66     67592 $name ||= $opts->{name};
412 252   100     668 $N //= $opts->{N};
413 252   100     1177 $periodic ||= $opts->{periodic};
414 252   100     935 $params //= $opts->{params};
415 252 100 100     972 $params = [$params] if defined $params && !ref $params;
416              
417 252         696 $name =~ s/_per$//;
418              
419 252 100       641 my $windows = $periodic ? \%periodic_windows : \%symmetric_windows;
420 252 100       877 unless ( $windows->{$name}) {
421 1 50       3 my $type = $periodic ? 'periodic' : 'symmetric';
422 1         38 PDL::Core::barf "window: Unknown $type window '$name'.";
423             }
424              
425 251         544 $self->{name} = $name;
426 251   100     656 $self->{N} = $N // die "Can't continue with undefined value for N";
427 250 100       624 die "Can't continue with zero value for N" if !$self->{N};
428 249         415 $self->{periodic} = $periodic;
429 249         514 $self->{params} = $params;
430 249 100       2347 $self->{code} = __PACKAGE__->can( $name . ( $periodic ? '_per' : '' ) );
431 249         1264 $self->{samples} = undef;
432 249         469 $self->{modfreqs} = undef;
433              
434 249         891 return $self;
435             }
436              
437             =head2 samples
438              
439             =for usage
440              
441             $win->samples;
442              
443             =for ref
444              
445             Generate and return a reference to the ndarray of C<$N> samples for the window
446             C<$win>. This is the real-space representation of the window.
447              
448             The samples are stored in the object C<$win>, but are regenerated every time
449             C is invoked. See the method L below.
450              
451             =for example
452              
453             For example:
454              
455             my $win = PDL::DSP::Windows->new( 8, 'hann' );
456             print $win->samples, "\n";
457              
458             =cut
459              
460             sub samples {
461 239     239 1 396 my $self = shift;
462 239   100     412 my @args = ( $self->{N}, @{ $self->{params} // [] } );
  239         971  
463 239         823 $self->{samples} = $self->{code}->(@args);
464             }
465              
466             =head2 modfreqs
467              
468             =for usage
469              
470             $win->modfreqs;
471              
472             =for ref
473              
474             Generate and return a reference to the ndarray of the modulus of the fourier
475             transform of the samples for the window C<$win>.
476              
477             These values are stored in the object C<$win>, but are regenerated every time
478             C is invoked. See the method L below.
479              
480             =head3 options
481              
482             =over
483              
484             =item min_bins => MIN
485              
486             This sets the minimum number of frequency bins. Defaults to 1000. If necessary,
487             the ndarray of window samples are padded with zeroes before the fourier transform
488             is performed.
489              
490             =back
491              
492             =cut
493              
494             sub modfreqs {
495 14     14 1 2029 require PDL::FFT;
496 14         23199 my $self = shift;
497 14         104 my %opts = PDL::Options::iparse( { min_bins => 1000 }, PDL::Options::ifhref(shift) );
498              
499 14         4657 my $data = $self->get_samples;
500              
501 14         2850 my $n = $data->nelem;
502 14 50       68 my $fn = $n > $opts{min_bins} ? 2 * $n : $opts{min_bins};
503              
504 14         32 $n--;
505              
506 14         58 my $freq = PDL::Core::zeroes($fn);
507 14         464 $freq->slice("0:$n") .= $data;
508              
509 14         863 PDL::FFT::realfft($freq);
510              
511 14         5591 my $real = PDL::Core::zeroes($freq);
512 14         41135 my $img = PDL::Core::zeroes($freq);
513 14         36263 my $mid = ( $freq->nelem ) / 2 - 1;
514 14         72 my $mid1 = $mid + 1;
515              
516 14         121 $real->slice("0:$mid") .= $freq->slice("$mid:0:-1");
517 14         1263 $real->slice("$mid1:-1") .= $freq->slice("0:$mid");
518 14         591 $img->slice("0:$mid") .= $freq->slice("-1:$mid1:-1");
519 14         679 $img->slice("$mid1:-1") .= $freq->slice("$mid1:-1");
520              
521 14         605 return $self->{modfreqs} = $real ** 2 + $img ** 2;
522             }
523              
524             =head2 get
525              
526             =for usage
527              
528             my $windata = $win->get('samples');
529              
530             =for ref
531              
532             Get an attribute (or list of attributes) of the window C<$win>. If attribute
533             C is requested, then the samples are created with the method
534             L if they don't exist.
535              
536             =for example
537              
538             For example:
539              
540             my $win = PDL::DSP::Windows->new( 8, 'hann' );
541             print $win->get('samples'), "\n";
542              
543             =cut
544              
545             sub get {
546 5     5 1 2333 my $self = shift;
547 5         11 my @res;
548              
549 5         19 for (@_) {
550 7 100       32 if ($_ eq 'samples') {
    100          
551 2         8 push @res, $self->get_samples;
552             }
553             elsif ($_ eq 'modfreqs') {
554 2         10 push @res, $self->get_modfreqs;
555             }
556             else {
557 3         8 push @res, $self->{$_};
558             }
559             };
560              
561 5 100       620 return wantarray ? @res : $res[0];
562             }
563              
564             =head2 get_samples
565              
566             =for usage
567              
568             my $windata = $win->get_samples;
569              
570             =for ref
571              
572             Return a reference to the pdl of samples for the Window instance C<$win>. The
573             samples will be generated with the method L if and only if they have
574             not yet been generated.
575              
576             =cut
577              
578             sub get_samples {
579 47     47 1 2534 my $self = shift;
580 47 100       211 return $self->{samples} if defined $self->{samples};
581 29         112 return $self->samples;
582             }
583              
584             =head2 get_modfreqs
585              
586             =for usage
587              
588             my $winfreqs = $win->get_modfreqs;
589             my $winfreqs = $win->get_modfreqs(OPTS);
590              
591             =for ref
592              
593             Return a reference to the pdl of the frequency response (modulus of the DFT)
594             for the Window instance C<$win>.
595              
596             Options passed as a hash reference will be passed to the L. The
597             data are created with L if they don't exist. The data are also
598             created even if they already exist if options are supplied. Otherwise the
599             cached data are returned.
600              
601             =head3 options
602              
603             =over
604              
605             =item min_bins => MIN
606              
607             This sets the minimum number of frequency bins. See L. Defaults
608             to 1000.
609              
610             =back
611              
612             =cut
613              
614             sub get_modfreqs {
615 13     13 1 4235 my $self = shift;
616 13 100       70 return $self->modfreqs(@_) if @_;
617 5 100       28 return $self->{modfreqs} if defined $self->{modfreqs};
618 2         9 return $self->modfreqs;
619             }
620              
621             =head2 get_params
622              
623             =for usage
624              
625             my $params = $win->get_params;
626              
627             =for ref
628              
629             Create a new array containing the parameter values for the instance C<$win>
630             and return a reference to the array. Note that not all window types take
631             parameters.
632              
633             =cut
634              
635 2     2 1 12 sub get_params { shift->{params} }
636              
637 1     1 0 9 sub get_N { shift->{N} }
638              
639             =head2 get_name
640              
641             =for usage
642              
643             print $win->get_name, "\n";
644              
645             =for ref
646              
647             Return a name suitable for printing associated with the window C<$win>. This is
648             something like the name used in the documentation for the particular window
649             function. This is static data and does not depend on the instance.
650              
651             =cut
652              
653             sub get_name {
654 29     29 1 11353 my $self = shift;
655              
656 29 100       184 if ( my $name = $window_print_names{ $self->{name} } ) {
657 1         7 return "$name window";
658             }
659              
660 28 100       110 if ( my $name = $window_names{ $self->{name} } ) {
661 6 100       36 return "$name window" unless $name =~ /^\*/;
662 5         57 return $name;
663             }
664              
665 22         156 return ucfirst $self->{name} . ' window';
666             }
667              
668             sub get_param_names {
669 6     6 0 14 my $self = shift;
670 6         22 my $params = $window_parameters{ $self->{name} };
671 6 100       42 return undef unless $params;
672 5         11 return [ @{$params} ];
  5         36  
673             }
674              
675             sub format_param_vals {
676 10     10 0 26 my $self = shift;
677              
678 10 100       22 my @params = @{ $self->{params} || [] };
  10         60  
679 10 100       68 return '' unless @params;
680              
681 4 50       12 my @names = @{ $self->get_param_names || [] };
  4         14  
682 4 50       16 return '' unless @names;
683              
684             join ', ', map {
685 4         16 my $name = $names[$_];
  10         22  
686 10         33 $name =~ s/^\$//;
687 10         59 join' = ', $name, $params[$_];
688             } 0 .. $#params;
689             }
690              
691             sub format_plot_param_vals {
692 8     8 0 26 my $ps = shift->format_param_vals;
693 8 100       55 return $ps ? ": $ps" : $ps;
694             }
695              
696             =head2 plot
697              
698             =for usage
699              
700             $win->plot;
701             $win->plot($pgswin); # can supply e.g. for multi-plotting
702              
703             =for ref
704              
705             Plot the samples. Uses L. The
706             default display type is used.
707              
708             =cut
709              
710             sub plot {
711 2 50   2 1 34 PDL::Core::barf 'PDL::DSP::Windows::plot PDL::Graphics::Simple not available!' unless eval { require PDL::Graphics::Simple };
  2         22  
712 2         4 my $self = shift;
713 2 50       10 my $pgsw = UNIVERSAL::isa($_[0], 'PDL::Graphics::Simple') ? shift : undef;
714 2         9 my @args = (
715             with => 'lines',
716             $self->get_samples,
717             { title => $self->get_name . $self->format_plot_param_vals,
718             xlabel => 'Time (samples)',
719             ylabel => 'amplitude' },
720             );
721 2 50       14 $pgsw ? $pgsw->plot(@args) : PDL::Graphics::Simple::plot(@args);
722 2         31 return $self;
723             }
724              
725             =head2 plot_freq
726              
727             =for usage
728              
729             Can be called like this
730              
731             $win->plot_freq;
732              
733             Or this
734              
735             $win->plot_freq({ ordinate => ORDINATE });
736             $win->plot_freq($pgswin, { ordinate => ORDINATE }); # can supply e.g. for multi-plotting
737              
738             =for ref
739              
740             Plot the frequency response (magnitude of the DFT of the window samples).
741             The response is plotted in dB, and the frequency (by default) as a fraction of
742             the Nyquist frequency. Uses L.
743             The default display type is used.
744              
745             =head3 options
746              
747             =over
748              
749             =item coord => COORD
750              
751             This sets the units of frequency of the co-ordinate axis. C must be one
752             of C, for fraction of the nyquist frequency (range C<-1, 1>);
753             C, for fraction of the sampling frequncy (range C<-0.5, 0.5>); or
754             C for frequency bin number (range C<0, $N - 1>). The default value is
755             C.
756              
757             =item min_bins => MIN
758              
759             This sets the minimum number of frequency bins. See L.
760             Defaults to 1000.
761              
762             =back
763              
764             =cut
765              
766             my %coord2xlab = (
767             nyquist=>'Fraction of Nyquist frequency',
768             sample=>'Fraction of sampling frequency',
769             bin=>'bin',
770             );
771             sub plot_freq {
772 6 50   6 1 64 PDL::Core::barf 'PDL::DSP::Windows::plot PDL::Graphics::Simple not available!' unless eval { require PDL::Graphics::Simple };
  6         82  
773 6         55 my $self = shift;
774 6 50       31 my $pgsw = UNIVERSAL::isa($_[0], 'PDL::Graphics::Simple') ? shift : undef;
775 6 50       38 my $opts = new PDL::Options({
776             coord => 'nyquist',
777             min_bins => 1000
778             })->options( @_ ? shift : {} );
779 6         1730 my $mf = $self->get_modfreqs({ min_bins => $opts->{min_bins} });
780 6         1407 $mf /= $mf->max;
781 6         455 my $title = $self->get_name . $self->format_plot_param_vals
782             . ', frequency response. ENBW=' . sprintf( '%2.3f', $self->enbw );
783 6         1139 my $coord = $opts->{coord};
784 6   66     37 my $xlab = $coord2xlab{$coord} //
785             PDL::Core::barf "plot_freq: Unknown ordinate unit specification $coord";
786 5         12 my $ylab = 'frequency response (dB)';
787             my $coordinate_range = $coord eq 'nyquist' ? 1 :
788             $coord eq 'sample' ? 0.5 :
789 5 100       29 $self->{N} / 2;
    100          
790 5         18 my $coordinates = PDL::Core::zeroes($mf)
791             ->xlinvals( -$coordinate_range, $coordinate_range );
792 5         13168 my @args = (
793             with => 'lines',
794             $coordinates,
795             20 * PDL::Ops::log10($mf),
796             { title => $title,
797             xrange => [-$coordinate_range,$coordinate_range],
798             xlabel => $xlab,
799             ylabel => $ylab },
800             );
801 5 50       628 $pgsw ? $pgsw->plot(@args) : PDL::Graphics::Simple::plot(@args);
802 5         118 return $self;
803             }
804              
805             =head2 enbw
806              
807             =for usage
808              
809             $win->enbw;
810              
811             =for ref
812              
813             Compute and return the equivalent noise bandwidth of the window.
814              
815             =cut
816              
817             sub enbw {
818 23     23 1 326 my $w = shift->get_samples;
819 23         27810 $w->nelem * ( $w ** 2 )->sum / $w->sum ** 2;
820             }
821              
822             =head2 coherent_gain
823              
824             =for usage
825              
826             $win->coherent_gain;
827              
828             =for ref
829              
830             Compute and return the coherent gain (the dc gain) of the window. This is just
831             the average of the samples.
832              
833             =cut
834              
835             sub coherent_gain {
836 1     1 1 11 my $w = shift->get_samples;
837 1         473 $w->sum / $w->nelem;
838             }
839              
840              
841             =head2 process_gain
842              
843             =for usage
844              
845             $win->coherent_gain;
846              
847             =for ref
848              
849             Compute and return the processing gain (the dc gain) of the window. This is
850             just the multiplicative inverse of the C.
851              
852             =cut
853              
854 1     1 1 1368 sub process_gain { 1 / shift->enbw }
855              
856             # not quite correct for some reason.
857             # Seems like 10*log10(this) / 1.154
858             # gives the correct answer in decibels
859              
860             =head2 scallop_loss
861              
862             =for usage
863              
864             $win->scallop_loss;
865              
866             =for ref
867              
868             Compute and return the scalloping loss of the window.
869              
870             =cut
871              
872             sub scallop_loss {
873 37     37 1 95 my $w = shift->samples;
874              
875             # Adapted from https://stackoverflow.com/a/40912607
876 37         8763 my $num = $w * exp( -( I() * PDL::sequence($w) * PI / $w->nelem ) );
877              
878 37         12374 20 * PDL::Ops::log10( abs( $num->sum ) / abs($w)->sum );
879             }
880              
881             =head1 WINDOW FUNCTIONS
882              
883             These window-function subroutines return a pdl of C<$N> samples. For most
884             windows, there are a symmetric and a periodic version. The symmetric versions
885             are functions of C<$N> points, uniformly spaced, and taking values from x_lo
886             through x_hi. Here, a periodic function of C< $N > points is equivalent to its
887             symmetric counterpart of C<$N + 1> points, with the final point omitted. The
888             name of a periodic window-function subroutine is the same as that for the
889             corresponding symmetric function, except it has the suffix C<_per>. The
890             descriptions below describe the symmetric version of each window.
891              
892             The term 'Blackman-Harris family' is meant to include the Hamming family and
893             the Blackman family. These are functions of sums of cosines.
894              
895             Unless otherwise noted, the arguments in the cosines of all symmetric window
896             functions are multiples of C<$N> numbers uniformly spaced from 0 through
897             2Ï€.
898              
899             =cut
900              
901             sub bartlett {
902 4 100   4 1 1630 PDL::Core::barf 'bartlett: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
903 2         7 my ($N) = @_;
904 2         10 1 - abs( PDL::Core::zeroes($N)->xlinvals( -1, 1 ) );
905             }
906              
907             sub bartlett_per {
908 4 100   4 0 1772 PDL::Core::barf 'bartlett: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
909 2         5 my ($N) = @_;
910 2         8 1 - abs( PDL::Core::zeroes($N)->xlinvals( -1, ( -1 + 1 * ( $N - 1 ) ) / $N ) );
911             }
912              
913             sub bartlett_hann {
914 6 100   6 1 2223 PDL::Core::barf 'bartlett_hann: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
915 4         11 my ($N) = @_;
916 4         20 0.62 - 0.48 * abs( PDL::Core::zeroes($N)->xlinvals( -0.5, 0.5 ) )
917             + 0.38 * cos( PDL::Core::zeroes($N)->xlinvals( - PI, PI ) );
918             }
919              
920             sub bartlett_hann_per {
921 4 100   4 0 2089 PDL::Core::barf 'bartlett_hann: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
922 2         18 my ($N) = @_;
923 2         10 0.62 - 0.48 * abs( PDL::Core::zeroes($N)->xlinvals( -0.5, ( -0.5 + 0.5 * ( $N - 1 ) ) / $N ) )
924             + 0.38 * cos( PDL::Core::zeroes($N)->xlinvals( - PI, ( - PI + PI * ( $N - 1 ) ) / $N ) );
925             }
926              
927             sub blackman {
928 6 100   6 1 1787 PDL::Core::barf 'blackman: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
929 4         11 my ($N) = @_;
930 4         21 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI ) );
931 4         1894 0.34 + $cx * ( -0.5 + $cx * 0.16 );
932             }
933              
934             sub blackman_per {
935 4 100   4 0 1886 PDL::Core::barf 'blackman: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
936 2         6 my ($N) = @_;
937 2         8 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI * ( $N - 1 ) / $N ) );
938 2         384 0.34 + $cx * ( -0.5 + $cx * 0.16 );
939             }
940              
941             sub blackman_bnh {
942 4 100   4 1 1957 PDL::Core::barf 'blackman_bnh: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
943 2         5 my ($N) = @_;
944 2         10 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI ) );
945 2         386 0.3461008 + $cx * ( -0.4973406 + $cx * 0.1565586 );
946             }
947              
948             sub blackman_bnh_per {
949 4 100   4 0 1836 PDL::Core::barf 'blackman_bnh: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
950 2         5 my ($N) = @_;
951 2         8 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI * ( $N - 1 ) / $N ) );
952 2         348 0.3461008 + $cx * ( -0.4973406 + $cx * 0.1565586 );
953             }
954              
955             sub blackman_ex {
956 5 100   5 1 1748 PDL::Core::barf 'blackman_ex: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
957 3         9 my ($N) = @_;
958 3         13 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI ) );
959 3         717 0.349742046431642 + $cx * ( -0.496560619088564 + $cx * 0.153697334479794 );
960             }
961              
962             sub blackman_ex_per {
963 4 100   4 0 1701 PDL::Core::barf 'blackman_ex: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
964 2         6 my ($N) = @_;
965 2         9 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI * ( $N - 1 ) / $N ) );
966 2         408 0.349742046431642 + $cx * ( -0.496560619088564 + $cx * 0.153697334479794 );
967             }
968              
969             sub blackman_gen {
970 8 100   8 1 3239 PDL::Core::barf 'blackman_gen: 2 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 2;
971 5         16 my ( $N, $alpha ) = @_;
972 5         24 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI ) );
973 5         1201 0.5 - $alpha + $cx * ( -0.5 + $cx * $alpha );
974             }
975              
976             sub blackman_gen_per {
977 5 100   5 0 2562 PDL::Core::barf 'blackman_gen: 2 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 2;
978 2         6 my ($N,$alpha) = @_;
979 2         9 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI * ( $N - 1 ) / $N ) );
980 2         348 0.5 - $alpha + $cx * ( -0.5 + $cx * $alpha );
981             }
982              
983             sub blackman_gen3 {
984 9 100   9 1 5419 PDL::Core::barf 'blackman_gen3: 4 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 4;
985 4         12 my ($N,$a0,$a1,$a2) = @_;
986 4         23 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI ) );
987 4         997 $a0 - $a2 + ( $cx * ( -$a1 + $cx * 2 * $a2 ) );
988             }
989              
990             sub blackman_gen3_per {
991 7 100   7 0 4622 PDL::Core::barf 'blackman_gen3: 4 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 4;
992 2         7 my ( $N, $a0, $a1, $a2 ) = @_;
993 2         7 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI * ( $N - 1 ) / $N ) );
994 2         342 $a0 - $a2 + ( $cx * ( -$a1 + $cx * 2 * $a2 ) );
995             }
996              
997             sub blackman_gen4 {
998 8 100   8 1 5390 PDL::Core::barf 'blackman_gen4: 5 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 5;
999 2         6 my ( $N, $a0, $a1, $a2, $a3 ) = @_;
1000 2         49 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI ) );
1001 2         388 $a0 - $a2 + $cx * ( ( -$a1 + 3 * $a3 ) + $cx * ( 2 * $a2 + $cx * -4 * $a3 ) );
1002             }
1003              
1004             sub blackman_gen4_per {
1005 8 100   8 0 5368 PDL::Core::barf 'blackman_gen4: 5 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 5;
1006 2         7 my ( $N, $a0, $a1, $a2, $a3 ) = @_;
1007 2         10 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI * ( $N - 1 ) / $N ) );
1008 2         396 $a0 - $a2 + $cx * ( ( -$a1 + 3 * $a3 ) + $cx * ( 2 * $a2 + $cx * -4 * $a3 ) );
1009             }
1010              
1011             sub blackman_gen5 {
1012 9 100   9 1 6651 PDL::Core::barf 'blackman_gen5: 6 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 6;
1013 2         7 my ( $N, $a0, $a1, $a2, $a3, $a4 ) = @_;
1014 2         8 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI ) );
1015 2         466 $a0 - $a2 + $a4 + $cx * ( ( -$a1 + 3 * $a3 )
1016             + $cx * ( 2 * $a2 - 8 * $a4 + $cx * ( -4 * $a3 + $cx * 8 * $a4 ) ) );
1017             }
1018              
1019             sub blackman_gen5_per {
1020 9 100   9 0 8660 PDL::Core::barf 'blackman_gen5: 6 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 6;
1021 2         8 my ( $N, $a0, $a1, $a2, $a3, $a4 ) = @_;
1022 2         8 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI * ( $N - 1 ) / $N ) );
1023 2         399 $a0 - $a2 + $a4 + $cx * ( ( -$a1 + 3 * $a3 )
1024             + $cx * ( 2 * $a2 - 8 * $a4 + $cx * ( -4 * $a3 + $cx * 8 * $a4 ) ) );
1025             }
1026              
1027             sub blackman_harris {
1028 5 100   5 1 1701 PDL::Core::barf 'blackman_harris: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1029 3         8 my ($N) = @_;
1030 3         14 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI ) );
1031 3         688 0.343103 + $cx * ( -0.49755 + $cx * 0.15844 );
1032             }
1033              
1034             sub blackman_harris_per {
1035 4 100   4 0 2334 PDL::Core::barf 'blackman_harris: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1036 2         6 my ($N) = @_;
1037 2         9 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI * ( $N - 1 ) / $N ) );
1038 2         538 0.343103 + $cx * ( -0.49755 + $cx * 0.15844 );
1039             }
1040              
1041             sub blackman_harris4 {
1042 6 100   6 1 1683 PDL::Core::barf 'blackman_harris4: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1043 4         11 my ($N) = @_;
1044 4         18 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI ) );
1045 4         1632 0.21747 + $cx * ( -0.45325 + $cx * ( 0.28256 + $cx * -0.04672 ) );
1046             }
1047              
1048             sub blackman_harris4_per {
1049 4 100   4 0 2013 PDL::Core::barf 'blackman_harris4: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1050 2         6 my ($N) = @_;
1051 2         9 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI * ( $N - 1 ) / $N ) );
1052 2         406 0.21747 + $cx * ( -0.45325 + $cx * ( 0.28256 + $cx * -0.04672 ) );
1053             }
1054              
1055             sub blackman_nuttall {
1056 6 100   6 1 1662 PDL::Core::barf 'blackman_nuttall: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1057 4         12 my ($N) = @_;
1058 4         17 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI ) );
1059 4         827 0.2269824 + $cx * ( -0.4572542 + $cx * ( 0.273199 + $cx * -0.0425644 ) );
1060             }
1061              
1062             sub blackman_nuttall_per {
1063 4 100   4 0 1655 PDL::Core::barf 'blackman_nuttall: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1064 2         4 my ($N) = @_;
1065 2         10 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI * ( $N - 1 ) / $N ) );
1066 2         392 0.2269824 + $cx * ( -0.4572542 + $cx * ( 0.273199 + $cx * -0.0425644 ) );
1067             }
1068              
1069             sub bohman {
1070 7 100   7 1 1815 PDL::Core::barf 'bohman: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1071 5         14 my ($N) = @_;
1072 5         22 my $x = abs( PDL::Core::zeroes($N)->xlinvals( -1, 1 ) );
1073 5         1409 ( 1 - $x ) * cos( PI * $x ) + ( 1 / PI ) * sin( PI * $x );
1074             }
1075              
1076             sub bohman_per {
1077 4 100   4 0 1872 PDL::Core::barf 'bohman: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1078 2         6 my ($N) = @_;
1079 2         7 my $x = abs( PDL::Core::zeroes($N)->xlinvals( -1, ( -1 + 1 * ( $N - 1 ) ) / $N ) );
1080 2         411 ( 1 - $x ) * cos( PI * $x ) + ( 1 / PI ) * sin( PI * $x );
1081             }
1082              
1083             sub cauchy {
1084 6 100   6 1 2662 PDL::Core::barf 'cauchy: 2 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 2;
1085 3         10 my ( $N, $alpha ) = @_;
1086 3         13 1 / ( 1 + ( PDL::Core::zeroes($N)->xlinvals( -1, 1 ) * $alpha ) ** 2 );
1087             }
1088              
1089             sub cauchy_per {
1090 5 100   5 0 2791 PDL::Core::barf 'cauchy: 2 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 2;
1091 2         6 my ( $N, $alpha ) = @_;
1092 2         9 1 / ( 1 + ( PDL::Core::zeroes($N)->xlinvals( -1, ( -1 + 1 * ( $N - 1 ) ) / $N ) * $alpha ) ** 2 );
1093             }
1094              
1095             sub chebyshev {
1096 8     8 1 3253 require PDL::FFT;
1097 8 100       11118 PDL::Core::barf 'chebyshev: 2 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 2;
1098 5         16 my ( $N, $at ) = @_;
1099              
1100 5         15 my ( $M, $M1, $pos, $pos1 );
1101              
1102 5         284 my $beta = PDL::Math::cosh( 1 / ( $N - 1 ) * PDL::Math::acosh( 1 / ( 10 ** ( -$at / 20 ) ) ) );
1103 5         257 my $x = $beta * cos( PI * PDL::Basic::sequence($N) / $N );
1104              
1105 5         916 my $cw = chebpoly( $N - 1, $x );
1106              
1107 5 100       51 if ( $N % 2 ) { # odd
1108 1         4 $M1 = ( $N + 1 ) / 2;
1109 1         3 $M = $M1 - 1;
1110 1         3 $pos = 0;
1111 1         2 $pos1 = 1;
1112              
1113 1         23 PDL::FFT::realfft($cw);
1114             }
1115             else { # half-sample delay (even order)
1116 4         22 my $arg = PI / $N * PDL::Basic::sequence($N);
1117 4         445 my $cw_im = $cw * sin($arg);
1118 4         259 $cw *= cos($arg);
1119              
1120 4         268 PDL::FFT::ifftnd( $cw, $cw_im );
1121              
1122 4         1304 $M1 = $N / 2;
1123 4         12 $M = $M1 - 1;
1124 4         10 $pos = 1;
1125 4         20 $pos1 = 0;
1126             }
1127              
1128 5         249 $cw /= $cw->at($pos);
1129              
1130 5         236 my $cwout = PDL::Core::zeroes($N);
1131              
1132 5         89 $cwout->slice("0:$M") .= $cw->slice("$M:0:-1");
1133 5         303 $cwout->slice("$M1:-1") .= $cw->slice("$pos1:$M");
1134 5         232 $cwout /= PDL::Ufunc::max($cwout);
1135              
1136 5         238 $cwout;
1137             }
1138              
1139             sub cos_alpha {
1140 10 100   10 1 2524 PDL::Core::barf 'cos_alpha: 2 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 2;
1141 7         21 my ( $N, $alpha ) = @_;
1142 7         32 ( sin( PDL::Core::zeroes($N)->xlinvals( 0, PI ) ) ) ** $alpha ;
1143             }
1144              
1145             sub cos_alpha_per {
1146 6 100   6 0 2779 PDL::Core::barf 'cos_alpha: 2 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 2;
1147 3         8 my ( $N, $alpha ) = @_;
1148 3         15 sin( PDL::Core::zeroes($N)->xlinvals( 0, PI * ( $N - 1 ) / $N ) ) ** $alpha ;
1149             }
1150              
1151             sub cosine {
1152 6 100   6 1 1839 PDL::Core::barf 'cosine: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1153 4         13 my ($N) = @_;
1154 4         35 sin( PDL::Core::zeroes($N)->xlinvals( 0, PI ) );
1155             }
1156              
1157             sub cosine_per {
1158 4 100   4 0 1834 PDL::Core::barf 'cosine: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1159 2         7 my ($N) = @_;
1160 2         9 sin( PDL::Core::zeroes($N)->xlinvals( 0, PI * ( $N - 1 ) / $N ) );
1161             }
1162              
1163             sub dpss {
1164 0 0   0 1 0 PDL::Core::barf 'dpss: 2 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 2;
1165 0         0 my ( $N, $beta ) = @_;
1166              
1167 0 0       0 PDL::Core::barf 'dpss: PDL::LinearAlgebra not installed.' unless eval { require PDL::LinearAlgebra::Special };
  0         0  
1168 0 0 0     0 PDL::Core::barf "dpss: $beta not between 0 and $N." unless $beta >= 0 and $beta <= $N;
1169              
1170 0         0 $beta /= $N / 2;
1171              
1172 0         0 my $k = PDL::Basic::sequence($N);
1173 0         0 my $s = sin( PI * $beta * $k ) / $k;
1174              
1175 0         0 $s->slice('0') .= $beta;
1176              
1177 0         0 my ( $ev, $e ) = PDL::MatrixOps::eigens_sym( PDL::LinearAlgebra::Special::mtoeplitz($s) );
1178 0         0 my $i = $e->maximum_ind;
1179              
1180 0         0 $ev->slice("($i)")->copy;
1181             }
1182              
1183             sub dpss_per {
1184 0 0   0 0 0 PDL::Core::barf 'dpss: 2 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 2;
1185 0         0 my ( $N, $beta ) = @_;
1186 0         0 $N++;
1187              
1188 0 0       0 PDL::Core::barf 'dpss: PDL::LinearAlgebra not installed.' unless eval { require PDL::LinearAlgebra::Special };
  0         0  
1189 0 0 0     0 PDL::Core::barf "dpss: $beta not between 0 and $N." unless $beta >= 0 and $beta <= $N;
1190              
1191 0         0 $beta /= $N / 2;
1192              
1193 0         0 my $k = PDL::Basic::sequence($N);
1194 0         0 my $s = sin( PI * $beta * $k ) / $k;
1195              
1196 0         0 $s->slice('0') .= $beta;
1197              
1198 0         0 my ( $ev, $e ) = PDL::MatrixOps::eigens_sym( PDL::LinearAlgebra::Special::mtoeplitz($s) );
1199 0         0 my $i = $e->maximum_ind;
1200              
1201 0         0 $ev->slice("($i),0:-2")->copy;
1202             }
1203              
1204             sub exponential {
1205 4 100   4 1 1869 PDL::Core::barf 'exponential: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1206 2         7 my ($N) = @_;
1207 2         8 2 ** ( 1 - abs( PDL::Core::zeroes($N)->xlinvals( -1, 1 ) ) ) - 1;
1208             }
1209              
1210             sub exponential_per {
1211 4 100   4 0 1894 PDL::Core::barf 'exponential: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1212 2         4 my ($N) = @_;
1213 2         9 2 ** ( 1 - abs( PDL::Core::zeroes($N)->xlinvals( -1, ( -1 + 1 * ( $N - 1 ) ) / $N ) ) ) - 1;
1214             }
1215              
1216             sub flattop {
1217 7 100   7 1 1857 PDL::Core::barf 'flattop: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1218 5         13 my ($N) = @_;
1219 5         23 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI ) );
1220 5         2005 -0.05473684 + $cx * ( -0.165894739 + $cx * ( 0.498947372 + $cx * ( -0.334315788 + $cx * 0.055578944 ) ) );
1221             }
1222              
1223             sub flattop_per {
1224 4 100   4 0 2123 PDL::Core::barf 'flattop: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1225 2         6 my ($N) = @_;
1226 2         9 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI * ( $N - 1 ) / $N ) );
1227 2         358 -0.05473684 + $cx * ( -0.165894739 + $cx * ( 0.498947372 + $cx * ( -0.334315788 + $cx * 0.055578944 ) ) );
1228             }
1229              
1230             sub gaussian {
1231 8 100   8 1 2945 PDL::Core::barf 'gaussian: 2 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 2;
1232 5         15 my ( $N, $beta ) = @_;
1233 5         25 exp( -0.5 * ( $beta * PDL::Core::zeroes($N)->xlinvals( -1, 1 ) ) ** 2);
1234             }
1235              
1236             sub gaussian_per {
1237 5 100   5 0 2475 PDL::Core::barf 'gaussian: 2 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 2;
1238 2         6 my ( $N, $beta ) = @_;
1239 2         10 exp( -0.5 * ( $beta * PDL::Core::zeroes($N)->xlinvals( -1, ( -1 + 1 * ( $N - 1 ) ) / $N ) ) ** 2 );
1240             }
1241              
1242             sub hamming {
1243 27 100   27 1 2062 PDL::Core::barf 'hamming: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1244 25         65 my ($N) = @_;
1245 25         132 0.54 + -0.46 * cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI ) );
1246             }
1247              
1248             sub hamming_per {
1249 4 100   4 0 1903 PDL::Core::barf 'hamming: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1250 2         4 my ($N) = @_;
1251 2         11 0.54 + -0.46 * cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI * ( $N - 1 ) / $N ) );
1252             }
1253              
1254             sub hamming_ex {
1255 4 100   4 1 1820 PDL::Core::barf 'hamming_ex: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1256 2         6 my ($N) = @_;
1257 2         10 0.53836 + -0.46164 * cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI ) );
1258             }
1259              
1260             sub hamming_ex_per {
1261 4 100   4 0 1815 PDL::Core::barf 'hamming_ex: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1262 2         5 my ($N) = @_;
1263 2         10 0.53836 + -0.46164 * cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI * ( $N - 1 ) / $N ) );
1264             }
1265              
1266             sub hamming_gen {
1267 5 100   5 1 3093 PDL::Core::barf 'hamming_gen: 2 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 2;
1268 2         6 my ( $N, $a ) = @_;
1269 2         10 $a - ( 1 - $a ) * cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI ) );
1270             }
1271              
1272             sub hamming_gen_per {
1273 5 100   5 0 2572 PDL::Core::barf 'hamming_gen: 2 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 2;
1274 2         6 my ( $N, $a ) = @_;
1275 2         18 $a - ( 1 - $a ) * cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI * ( $N - 1 ) / $N ) );
1276             }
1277              
1278             sub hann {
1279 8 100   8 1 1786 PDL::Core::barf 'hann: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1280 6         17 my ($N) = @_;
1281 6         30 0.5 + -0.5 * cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI ) );
1282             }
1283              
1284             sub hann_per {
1285 4 100   4 0 1822 PDL::Core::barf 'hann: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1286 2         6 my ($N) = @_;
1287 2         28 0.5 + -0.5 * cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI * ( $N - 1 ) / $N ) );
1288             }
1289              
1290             sub hann_matlab {
1291 3 100   3 1 1791 PDL::Core::barf 'hann_matlab: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1292 1         4 my ($N) = @_;
1293 1         6 0.5 - 0.5 * cos( PDL::Core::zeroes($N)->xlinvals( TPI / ( $N + 1 ), TPI * $N / ( $N + 1 ) ) );
1294             }
1295              
1296             sub hann_poisson {
1297 9 100   9 1 4505 PDL::Core::barf 'hann_poisson: 2 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 2;
1298 6         17 my ( $N, $alpha ) = @_;
1299 6         29 0.5 * ( 1 + cos( PDL::Core::zeroes($N)->xlinvals( - PI, PI ) ) )
1300             * exp( -$alpha * abs( PDL::Core::zeroes($N)->xlinvals( -1, 1 ) ) );
1301              
1302             }
1303              
1304             sub hann_poisson_per {
1305 5 100   5 0 2481 PDL::Core::barf 'hann_poisson: 2 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 2;
1306 2         6 my ( $N, $alpha ) = @_;
1307 2         48 0.5 * ( 1 + cos( PDL::Core::zeroes($N)->xlinvals( - PI, ( - PI + PI * ( $N - 1 ) ) / $N ) ) )
1308             * exp( -$alpha * abs( PDL::Core::zeroes($N)->xlinvals( -1, ( -1 + 1 * ( $N - 1 ) ) / $N ) ) );
1309             }
1310              
1311             sub kaiser {
1312 0 0   0 1 0 PDL::Core::barf 'kaiser: 2 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 2;
1313 0         0 my ( $N, $beta ) = @_;
1314              
1315 0 0       0 PDL::Core::barf 'kaiser: PDL::GSLSF not installed' unless eval { require PDL::GSLSF::BESSEL };
  0         0  
1316              
1317 0         0 $beta *= PI;
1318              
1319 0         0 my ($n) = PDL::GSLSF::BESSEL::gsl_sf_bessel_In(
1320             $beta * sqrt( 1 - PDL::Core::zeroes($N)->xlinvals( -1, 1 ) ** 2 ), 0 );
1321              
1322 0         0 my ($d) = PDL::GSLSF::BESSEL::gsl_sf_bessel_In( $beta, 0 );
1323              
1324 0         0 $n / $d;
1325             }
1326              
1327             sub kaiser_per {
1328 0 0   0 0 0 PDL::Core::barf 'kaiser: 2 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 2;
1329 0         0 my ($N,$beta) = @_;
1330              
1331 0 0       0 PDL::Core::barf 'kaiser: PDL::GSLSF not installed' unless eval { require PDL::GSLSF::BESSEL };
  0         0  
1332              
1333 0         0 $beta *= PI;
1334              
1335 0         0 my ($n) = PDL::GSLSF::BESSEL::gsl_sf_bessel_In(
1336             $beta * sqrt( 1 - PDL::Core::zeroes($N)->xlinvals( -1, ( -1 + 1 * ( $N - 1 ) ) / $N ) ** 2 ), 0);
1337              
1338 0         0 my ($d) = PDL::GSLSF::BESSEL::gsl_sf_bessel_In( $beta, 0 );
1339              
1340 0         0 $n / $d;
1341             }
1342              
1343             sub lanczos {
1344 6 100   6 1 1642 PDL::Core::barf 'lanczos: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1345 4         12 my ($N) = @_;
1346              
1347 4         18 my $x = PI * PDL::Core::zeroes($N)->xlinvals( -1, 1 );
1348 4         1086 my $res = sin($x) / $x;
1349              
1350 4 100       868 $res->slice( int( $N / 2 ) ) .= 1 if $N % 2;
1351              
1352 4         85 $res;
1353             }
1354              
1355             sub lanczos_per {
1356 4 100   4 0 1811 PDL::Core::barf 'lanczos: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1357 2         6 my ($N) = @_;
1358              
1359 2         9 my $x = PI * PDL::Core::zeroes($N)->xlinvals( -1, ( -1 + 1 * ( $N - 1 ) ) / $N );
1360 2         369 my $res = sin($x) / $x;
1361              
1362 2 100       66 $res->slice( int( $N / 2 ) ) .= 1 unless $N % 2;
1363              
1364 2         66 $res;
1365             }
1366              
1367             sub nuttall {
1368 5 100   5 1 1770 PDL::Core::barf 'nuttall: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1369 3         9 my ($N) = @_;
1370 3         16 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI ) );
1371 3         647 0.2269824 + $cx * ( -0.4572542 + $cx * ( 0.273199 + $cx * -0.0425644 ) );
1372             }
1373              
1374             sub nuttall_per {
1375 4 100   4 0 2011 PDL::Core::barf 'nuttall: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1376 2         7 my ($N) = @_;
1377 2         8 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI * ( $N - 1 ) / $N ) );
1378 2         384 0.2269824 + $cx * ( -0.4572542 + $cx * ( 0.273199 + $cx * -0.0425644 ) );
1379             }
1380              
1381             sub nuttall1 {
1382 4 100   4 1 1666 PDL::Core::barf 'nuttall1: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1383 2         6 my ($N) = @_;
1384 2         9 my $cx = (cos(PDL::Core::zeroes($N)->xlinvals(0,TPI)));
1385              
1386 2         439 (0.211536) + ($cx * ((-0.449584) + ($cx * (0.288464 + $cx * (-0.050416) ))));
1387             }
1388              
1389             sub nuttall1_per {
1390 4 100   4 0 1948 PDL::Core::barf 'nuttall1: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1391 2         7 my ($N) = @_;
1392 2         8 my $cx = cos( PDL::Core::zeroes($N)->xlinvals( 0, TPI * ( $N - 1 ) / $N ) );
1393 2         399 0.211536 + $cx * ( -0.449584 + $cx * ( 0.288464 + $cx * -0.050416 ) );
1394             }
1395              
1396             sub parzen {
1397 7 100   7 1 1969 PDL::Core::barf 'parzen: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1398 5         16 my ($N) = @_;
1399              
1400 5         23 my $x = PDL::Core::zeroes($N)->xlinvals( -1, 1 );
1401              
1402 5         1350 my $x1 = $x->where( $x <= -0.5 );
1403 5         1075 my $x2 = $x->where( ( $x < 0.5 ) & ( $x > -0.5 ) );
1404 5         1485 my $x3 = $x->where( $x >= 0.5 );
1405              
1406 5         835 $x1 .= 2 * ( 1 - abs($x1) ) ** 3;
1407 5         1047 $x2 .= 1 - 6 * $x2 ** 2 * ( 1 - abs($x2) );
1408 5         1587 $x3 .= $x1->slice('-1:0:-1');
1409              
1410 5         358 $x;
1411             }
1412              
1413             sub parzen_per {
1414 4 100   4 0 1697 PDL::Core::barf 'parzen: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1415 2         5 my ($N) = @_;
1416              
1417 2         10 my $x = PDL::Core::zeroes($N)->xlinvals( -1, ( -1 + ( $N - 1 ) ) / $N);
1418              
1419 2         372 my $x1 = $x->where( $x <= -0.5 );
1420 2         306 my $x2 = $x->where( ( $x < 0.5 ) & ( $x > -0.5 ) );
1421 2         326 my $x3 = $x->where( $x >= 0.5 );
1422              
1423 2         222 $x1 .= 2 * ( 1 - abs($x1)) ** 3;
1424 2         203 $x2 .= 1 - 6 * $x2 ** 2 * ( 1 - abs($x2) );
1425 2         214 $x3 .= $x1->slice('-1:1:-1');
1426              
1427 2         92 $x;
1428             }
1429              
1430             sub parzen_octave {
1431 4 100   4 1 2140 PDL::Core::barf 'parzen_octave: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1432 2         5 my ($N) = @_;
1433              
1434 2         9 my $r = ( $N - 1 ) / 2;
1435 2         4 my $N2 = $N / 2;
1436 2         6 my $r4 = $r / 2;
1437 2         16 my $n = PDL::Basic::sequence( 2 * $r + 1 ) - $r;
1438              
1439 2         365 my $n1 = $n->where( abs($n) <= $r4 );
1440 2         643 my $n2 = $n->where( $n > $r4 );
1441 2         482 my $n3 = $n->where( $n < -$r4 );
1442              
1443 2         452 $n1 .= 1 - 6 * ( abs($n1) / $N2 ) ** 2 + 6 * ( abs($n1) / $N2 ) ** 3;
1444 2         1926 $n2 .= 2 * ( 1 - abs($n2) / $N2 ) ** 3;
1445 2         607 $n3 .= 2 * ( 1 - abs($n3) / $N2 ) ** 3;
1446              
1447 2         631 $n;
1448             }
1449              
1450             sub poisson {
1451 9 100   9 1 3012 PDL::Core::barf 'poisson: 2 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 2;
1452 6         20 my ( $N, $alpha ) = @_;
1453 6         27 exp( -$alpha * abs( PDL::Core::zeroes($N)->xlinvals( -1, 1 ) ) );
1454             }
1455              
1456             sub poisson_per {
1457 5 100   5 0 2807 PDL::Core::barf 'poisson: 2 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 2;
1458 2         5 my ( $N, $alpha ) = @_;
1459 2         11 exp( -$alpha * abs( PDL::Core::zeroes($N)->xlinvals( -1, ( -1 + 1 * ( $N - 1 ) ) / $N ) ) );
1460             }
1461              
1462             sub rectangular {
1463 7 100   7 1 1774 PDL::Core::barf 'rectangular: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1464 5         15 my ($N) = @_;
1465 5         25 PDL::Core::ones($N);
1466             }
1467              
1468             sub rectangular_per {
1469 4 100   4 0 1851 PDL::Core::barf 'rectangular: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1470 2         33 my ($N) = @_;
1471 2         11 PDL::Core::ones($N);
1472             }
1473              
1474             sub triangular {
1475 7 100   7 1 1988 PDL::Core::barf 'triangular: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1476 5         15 my ($N) = @_;
1477 5         23 1 - abs( PDL::Core::zeroes($N)->xlinvals( -( $N - 1 ) / $N, ( $N - 1 ) / $N ) );
1478             }
1479              
1480             sub triangular_per {
1481 4 100   4 0 1676 PDL::Core::barf 'triangular: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1482 2         6 my ($N) = @_;
1483 2         10 1 - abs( PDL::Core::zeroes($N)->xlinvals( -$N / ( $N + 1 ), -1 / ( $N + 1 ) + ( $N - 1 ) / ( $N + 1 ) ) );
1484             }
1485              
1486             sub tukey {
1487 18 100   18 1 4305 PDL::Core::barf 'tukey: 2 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 2;
1488 15         29 my ( $N, $alpha ) = @_;
1489              
1490 15 100 100     75 PDL::Core::barf('tukey: alpha must be between 0 and 1') unless $alpha >= 0 and $alpha <= 1;
1491              
1492 13 50       31 return PDL::Core::ones($N) if $alpha == 0;
1493              
1494 13         50 my $x = PDL::Core::zeroes($N)->xlinvals( 0, 1 );
1495              
1496 13         2502 my $x1 = $x->where( $x < $alpha / 2 );
1497 13         2376 my $x2 = $x->where( ( $x <= 1 - $alpha / 2 ) & ( $x >= $alpha / 2 ) );
1498 13         2489 my $x3 = $x->where( $x > 1 - $alpha / 2 );
1499              
1500 13         1619 $x1 .= 0.5 * ( 1 + cos( PI * ( 2 * $x1 / $alpha - 1 ) ) );
1501 13         1677 $x2 .= 1;
1502 13         359 $x3 .= $x1->slice('-1:0:-1');
1503              
1504 13         561 return $x;
1505             }
1506              
1507             sub tukey_per {
1508 12 100   12 0 4184 PDL::Core::barf 'tukey: 2 arguments expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 2;
1509 9         19 my ( $N, $alpha ) = @_;
1510              
1511 9 100 100     50 PDL::Core::barf 'tukey: alpha must be between 0 and 1' unless $alpha >= 0 && $alpha <= 1;
1512              
1513 7 50       15 return PDL::Core::ones($N) if $alpha == 0;
1514              
1515 7         25 my $x = PDL::Core::zeroes($N)->xlinvals( 0, ( $N - 1 ) / $N );
1516              
1517 7         1088 my $x1 = $x->where( $x < $alpha / 2 );
1518 7         724 my $x2 = $x->where( ( $x <= 1 - $alpha / 2 ) & ( $x >= $alpha / 2 ) );
1519 7         809 my $x3 = $x->where( $x > 1 - $alpha / 2 );
1520              
1521 7         599 $x1 .= 0.5 * ( 1 + cos( PI * ( 2 * $x1 / $alpha - 1 ) ) );
1522 7         720 $x2 .= 1;
1523 7         123 $x3 .= $x1->slice('-1:1:-1');
1524              
1525 7         250 return $x;
1526             }
1527              
1528             sub welch {
1529 6 100   6 1 4016 PDL::Core::barf 'welch: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1530 4         10 my ($N) = @_;
1531 4         17 1 - PDL::Core::zeroes($N)->xlinvals( -1, 1 ) ** 2;
1532             }
1533              
1534             sub welch_per {
1535 4 100   4 0 1705 PDL::Core::barf 'welch: 1 argument expected. Got ' . scalar(@_) . ' arguments.' unless @_ == 1;
1536 2         6 my ($N) = @_;
1537 2         10 1 - PDL::Core::zeroes($N)->xlinvals( -1, ( -1 + 1 * ( $N - 1 ) ) / $N ) ** 2;
1538             }
1539              
1540             =head1 Symmetric window functions
1541              
1542             =head2 bartlett($N)
1543              
1544             The Bartlett window. (Ref 1). Another name for this window is the fejer window.
1545             This window is defined by
1546              
1547             1 - abs(arr)
1548              
1549             where the points in arr range from -1 through 1. See also
1550             L.
1551              
1552             =head2 bartlett_hann($N)
1553              
1554             The Bartlett-Hann window. Another name for this window is the Modified
1555             Bartlett-Hann window. This window is defined by
1556              
1557             0.62 - 0.48 * abs(arr) + 0.38 * arr1
1558              
1559             where the points in arr range from -1/2 through 1/2, and arr1 are the cos of
1560             points ranging from -PI through PI.
1561              
1562             =head2 blackman($N)
1563              
1564             The 'classic' Blackman window. (Ref 1). One of the Blackman-Harris family, with coefficients
1565              
1566             a0 = 0.42
1567             a1 = 0.5
1568             a2 = 0.08
1569              
1570             =head2 blackman_bnh($N)
1571              
1572             The Blackman-Harris (bnh) window. An improved version of the 3-term
1573             Blackman-Harris window given by Nuttall (Ref 2, p. 89). One of the
1574             Blackman-Harris family, with coefficients
1575              
1576             a0 = 0.4243801
1577             a1 = 0.4973406
1578             a2 = 0.0782793
1579              
1580             =head2 blackman_ex($N)
1581              
1582             The 'exact' Blackman window. (Ref 1). One of the Blackman-Harris family, with
1583             coefficients
1584              
1585             a0 = 0.426590713671539
1586             a1 = 0.496560619088564
1587             a2 = 0.0768486672398968
1588              
1589             =head2 blackman_gen($N,$alpha)
1590              
1591             The General classic Blackman window. A single parameter family of the 3-term
1592             Blackman window. This window is defined by
1593              
1594             my $cx = arr;
1595             .5 - $alpha + ($cx * (-.5 + $cx * $alpha));
1596              
1597             where the points in arr are the cos of points ranging from 0 through 2PI.
1598              
1599             =head2 blackman_gen3($N,$a0,$a1,$a2)
1600              
1601             The general form of the Blackman family. One of the Blackman-Harris family,
1602             with coefficients
1603              
1604             a0 = $a0
1605             a1 = $a1
1606             a2 = $a2
1607              
1608             =head2 blackman_gen4($N,$a0,$a1,$a2,$a3)
1609              
1610             The general 4-term Blackman-Harris window. One of the Blackman-Harris family,
1611             with coefficients
1612              
1613             a0 = $a0
1614             a1 = $a1
1615             a2 = $a2
1616             a3 = $a3
1617              
1618             =head2 blackman_gen5($N,$a0,$a1,$a2,$a3,$a4)
1619              
1620             The general 5-term Blackman-Harris window. One of the Blackman-Harris family,
1621             with coefficients
1622              
1623             a0 = $a0
1624             a1 = $a1
1625             a2 = $a2
1626             a3 = $a3
1627             a4 = $a4
1628              
1629             =head2 blackman_harris($N)
1630              
1631             The Blackman-Harris window. (Ref 1). One of the Blackman-Harris family, with
1632             coefficients
1633              
1634             a0 = 0.422323
1635             a1 = 0.49755
1636             a2 = 0.07922
1637              
1638             Another name for this window is the Minimum three term (sample) Blackman-Harris
1639             window.
1640              
1641             =head2 blackman_harris4($N)
1642              
1643             The minimum (sidelobe) four term Blackman-Harris window. (Ref 1). One of the
1644             Blackman-Harris family, with coefficients
1645              
1646             a0 = 0.35875
1647             a1 = 0.48829
1648             a2 = 0.14128
1649             a3 = 0.01168
1650              
1651             Another name for this window is the Blackman-Harris window.
1652              
1653             =head2 blackman_nuttall($N)
1654              
1655             The Blackman-Nuttall window. One of the Blackman-Harris family, with
1656             coefficients
1657              
1658             a0 = 0.3635819
1659             a1 = 0.4891775
1660             a2 = 0.1365995
1661             a3 = 0.0106411
1662              
1663             =head2 bohman($N)
1664              
1665             The Bohman window. (Ref 1). This window is defined by
1666              
1667             my $x = abs(arr);
1668             (1 - $x) * cos(PI * $x) + (1 / PI) * sin(PI * $x)
1669              
1670             where the points in arr range from -1 through 1.
1671              
1672             =head2 cauchy($N,$alpha)
1673              
1674             The Cauchy window. (Ref 1). Other names for this window are: Abel, Poisson.
1675             This window is defined by
1676              
1677             1 / (1 + (arr * $alpha) ** 2)
1678              
1679             where the points in arr range from -1 through 1.
1680              
1681             =head2 chebyshev($N,$at)
1682              
1683             The Chebyshev window. The frequency response of this window has C<$at> dB of
1684             attenuation in the stop-band. Another name for this window is the
1685             Dolph-Chebyshev window. No periodic version of this window is defined. This
1686             routine gives the same result as the routine B in Octave 3.6.2.
1687              
1688             =head2 cos_alpha($N,$alpha)
1689              
1690             The Cos_alpha window. (Ref 1). Another name for this window is the
1691             Power-of-cosine window. This window is defined by
1692              
1693             arr ** $alpha
1694              
1695             where the points in arr are the sin of points ranging from 0 through PI.
1696              
1697             =head2 cosine($N)
1698              
1699             The Cosine window. Another name for this window is the sine window. This
1700             window is defined by
1701              
1702             arr
1703              
1704             where the points in arr are the sin of points ranging from 0 through PI.
1705              
1706             =head2 dpss($N,$beta)
1707              
1708             The Digital Prolate Spheroidal Sequence (DPSS) window. The parameter C<$beta>
1709             is the half-width of the mainlobe, measured in frequency bins. This window
1710             maximizes the power in the mainlobe for given C<$N> and C<$beta>. Another
1711             name for this window is the sleppian window.
1712              
1713             =head2 exponential($N)
1714              
1715             The Exponential window. This window is defined by
1716              
1717             2 ** (1 - abs arr) - 1
1718              
1719             where the points in arr range from -1 through 1.
1720              
1721             =head2 flattop($N)
1722              
1723             The flat top window. One of the Blackman-Harris family, with coefficients
1724              
1725             a0 = 0.21557895
1726             a1 = 0.41663158
1727             a2 = 0.277263158
1728             a3 = 0.083578947
1729             a4 = 0.006947368
1730              
1731             =head2 gaussian($N,$beta)
1732              
1733             The Gaussian window. (Ref 1). Another name for this window is the Weierstrass
1734             window. This window is defined by
1735              
1736             exp (-0.5 * ($beta * arr )**2),
1737              
1738             where the points in arr range from -1 through 1.
1739              
1740             =head2 hamming($N)
1741              
1742             The Hamming window. (Ref 1). One of the Blackman-Harris family, with
1743             coefficients
1744              
1745             a0 = 0.54
1746             a1 = 0.46
1747              
1748             =head2 hamming_ex($N)
1749              
1750             The 'exact' Hamming window. (Ref 1). One of the Blackman-Harris family, with
1751             coefficients
1752              
1753             a0 = 0.53836
1754             a1 = 0.46164
1755              
1756             =head2 hamming_gen($N,$a)
1757              
1758             The general Hamming window. (Ref 1). One of the Blackman-Harris family, with
1759             coefficients
1760              
1761             a0 = $a
1762             a1 = (1-$a)
1763              
1764             =head2 hann($N)
1765              
1766             The Hann window. (Ref 1). One of the Blackman-Harris family, with coefficients
1767              
1768             a0 = 0.5
1769             a1 = 0.5
1770              
1771             Another name for this window is the hanning window. See also
1772             L.
1773              
1774             =head2 hann_matlab($N)
1775              
1776             The Hann (matlab) window. Equivalent to the Hann window of N+2 points, with the
1777             endpoints (which are both zero) removed. No periodic version of this window is
1778             defined. This window is defined by
1779              
1780             0.5 - 0.5 * arr
1781              
1782             where the points in arr are the cosine of points ranging from 2PI/($N+1)
1783             through 2PI*$N/($N+1). This routine gives the same result as the routine
1784             B in Matlab. See also L.
1785              
1786             =head2 hann_poisson($N,$alpha)
1787              
1788             The Hann-Poisson window. (Ref 1). This window is defined by
1789              
1790             0.5 * (1 + arr1) * exp (-$alpha * abs arr)
1791              
1792             where the points in arr range from -1 through 1, and arr1 are the cos of points
1793             ranging from -PI through PI.
1794              
1795             =head2 kaiser($N,$beta)
1796              
1797             The Kaiser window. (Ref 1). The parameter C<$beta> is the approximate
1798             half-width of the mainlobe, measured in frequency bins. Another name for this
1799             window is the Kaiser-Bessel window. This window is defined by
1800              
1801             $beta *= PI;
1802             my @n = gsl_sf_bessel_In($beta * sqrt(1 - arr ** 2), 0);
1803             my @d = gsl_sf_bessel_In($beta, 0);
1804             $n[0] / $d[0];
1805              
1806             where the points in arr range from -1 through 1.
1807              
1808             =head2 lanczos($N)
1809              
1810             The Lanczos window. Another name for this window is the sinc window. This
1811             window is defined by
1812              
1813             my $x = PI * arr;
1814             my $res = sin($x) / $x;
1815             my $mid;
1816             $mid = int($N / 2), $res->slice($mid) .= 1 if $N % 2;
1817             $res;
1818              
1819             where the points in arr range from -1 through 1.
1820              
1821             =head2 nuttall($N)
1822              
1823             The Nuttall window. One of the Blackman-Harris family, with coefficients
1824              
1825             a0 = 0.3635819
1826             a1 = 0.4891775
1827             a2 = 0.1365995
1828             a3 = 0.0106411
1829              
1830             See also L.
1831              
1832             =head2 nuttall1($N)
1833              
1834             The Nuttall (v1) window. A window referred to as the Nuttall window. One of the
1835             Blackman-Harris family, with coefficients
1836              
1837             a0 = 0.355768
1838             a1 = 0.487396
1839             a2 = 0.144232
1840             a3 = 0.012604
1841              
1842             This routine gives the same result as the routine B in Octave 3.6.2.
1843             See also L.
1844              
1845             =head2 parzen($N)
1846              
1847             The Parzen window. (Ref 1). Other names for this window are: Jackson,
1848             Valle-Poussin. This function disagrees with the Octave subroutine B,
1849             but agrees with Ref. 1. See also L.
1850              
1851             =head2 parzen_octave($N)
1852              
1853             The Parzen window. No periodic version of this window is defined. This routine
1854             gives the same result as the routine B in Octave 3.6.2. See also
1855             L.
1856              
1857             =head2 poisson($N,$alpha)
1858              
1859             The Poisson window. (Ref 1). This window is defined by
1860              
1861             exp(-$alpha * abs(arr))
1862              
1863             where the points in arr range from -1 through 1.
1864              
1865             =head2 rectangular($N)
1866              
1867             The Rectangular window. (Ref 1). Other names for this window are: dirichlet,
1868             boxcar.
1869              
1870             =head2 triangular($N)
1871              
1872             The Triangular window. This window is defined by
1873              
1874             1 - abs(arr)
1875              
1876             where the points in arr range from -$N/($N-1) through $N/($N-1).
1877             See also L.
1878              
1879             =head2 tukey($N,$alpha)
1880              
1881             The Tukey window. (Ref 1). Another name for this window is the tapered cosine
1882             window.
1883              
1884             =head2 welch($N)
1885              
1886             The Welch window. (Ref 1). Other names for this window are: Riez, Bochner,
1887             Parzen, parabolic. This window is defined by
1888              
1889             1 - arr**2
1890              
1891             where the points in arr range from -1 through 1.
1892              
1893             =head1 AUXILIARY SUBROUTINES
1894              
1895             These subroutines are used internally, but are also available for export.
1896              
1897             =head2 cos_mult_to_pow
1898              
1899             Convert Blackman-Harris coefficients. The BH windows are usually defined via
1900             coefficients for cosines of integer multiples of an argument. The same windows
1901             may be written instead as terms of powers of cosines of the same argument.
1902             These may be computed faster as they replace evaluation of cosines with
1903             multiplications. This subroutine is used internally to implement the
1904             Blackman-Harris family of windows more efficiently.
1905              
1906             This subroutine takes between 1 and 7 numeric arguments a0, a1, ...
1907              
1908             It converts the coefficients of this
1909              
1910             a0 - a1 cos(arg) + a2 cos(2 * arg) - a3 cos(3 * arg) + ...
1911              
1912             To the cofficients of this
1913              
1914             c0 + c1 cos(arg) + c2 cos(arg)**2 + c3 cos(arg)**3 + ...
1915              
1916             =head2 cos_pow_to_mult
1917              
1918             This function is the inverse of L.
1919              
1920             This subroutine takes between 1 and 7 numeric arguments c0, c1, ...
1921              
1922             It converts the coefficients of this
1923              
1924             c0 + c1 cos(arg) + c2 cos(arg)**2 + c3 cos(arg)**3 + ...
1925              
1926             To the cofficients of this
1927              
1928             a0 - a1 cos(arg) + a2 cos( 2 * arg) - a3 cos( 3 * arg) + ...
1929              
1930             =cut
1931              
1932             sub cos_pow_to_mult {
1933 9     9 1 19884 my @cin = @_;
1934 9 100       41 PDL::Core::barf 'cos_pow_to_mult: must be less than 8 args .' if @cin > 7;
1935 8         16 my $ex = 7 - @cin;
1936 8         27 my @c = ( @cin, (0) x $ex );
1937 8         38 my @as = (
1938             10 * $c[6] + 12 * $c[4] + 16 * $c[2] + 32 * $c[0],
1939             20 * $c[5] + 24 * $c[3] + 32 * $c[1],
1940             15 * $c[6] + 16 * $c[4] + 16 * $c[2],
1941             10 * $c[5] + 8 * $c[3],
1942             6 * $c[6] + 4 * $c[4],
1943             2 * $c[5],
1944             $c[6],
1945             );
1946 8 100       28 splice @as, -$ex if $ex;
1947 8         15 my $sign = -1;
1948 8         18 foreach (@as) {
1949 28         65 $_ /= -$sign * 32;
1950 28         55 $sign *= -1;
1951             }
1952 8         55 @as;
1953             }
1954              
1955             =head2 chebpoly
1956              
1957             =for usage
1958              
1959             chebpoly($n,$x)
1960              
1961             =for ref
1962              
1963             Returns the value of the C<$n>-th order Chebyshev polynomial at point C<$x>.
1964             $n and $x may be scalar numbers, pdl's, or array refs. However, at least one
1965             of $n and $x must be a scalar number.
1966              
1967             All mixtures of pdls and scalars could be handled much more easily as a PP
1968             routine. But, at this point PDL::DSP::Windows is pure perl/pdl, requiring no
1969             C/Fortran compiler.
1970              
1971             =cut
1972              
1973             sub chebpoly {
1974 8 50   8 1 461806 PDL::Core::barf 'chebpoly: Two arguments expected. Got ' .scalar(@_) . "\n" unless @_ == 2;
1975              
1976 8         29 my ( $n, $x ) = @_;
1977              
1978 8 100       50 if ( ref $x ) {
1979 7 50       23 PDL::Core::barf "chebpoly: neither $n nor $x is a scalar number" if ref($n);
1980 7         31 $x = PDL::Core::topdl($x);
1981              
1982 7         148 my $tn = PDL::Core::zeroes($x);
1983              
1984 7         7940 my ( $ind1, $ind2 ) = PDL::Primitive::which_both( abs($x) <= 1 );
1985              
1986 7         3253 $tn->index($ind1) .= cos( $n * PDL::Math::acos( $x->index($ind1) ) );
1987 7         1088 $tn->index($ind2) .= PDL::Math::cosh( $n * PDL::Math::acosh( $x->index($ind2) ) );
1988              
1989 7         415 return $tn;
1990             }
1991              
1992 1 50       4 $n = PDL::Core::topdl($n) if ref $n;
1993 1 50       35 return abs($x) <= 1 ? PDL::Math::cos( $n * PDL::Math::acos($x) ) : PDL::Math::cosh( $n * PDL::Math::acosh($x) );
1994             }
1995              
1996              
1997             sub cos_mult_to_pow {
1998 9     9 1 1671 my( @ain ) = @_;
1999 9 100       53 PDL::Core::barf 'cos_mult_to_pow: must be less than 8 args .' if @ain > 7;
2000 8         17 my $ex = 7 - @ain;
2001 8         21 my @a = ( @ain, (0) x $ex );
2002 8         47 my @cs = (
2003             -$a[6] + $a[4] - $a[2] + $a[0],
2004             -5 * $a[5] + 3 * $a[3] - $a[1],
2005             18 * $a[6] - 8 * $a[4] + 2 * $a[2],
2006             20 * $a[5] - 4 * $a[3],
2007             8 * $a[4] - 48 * $a[6],
2008             -16 * $a[5],
2009             32 * $a[6]
2010             );
2011 8 100       20 splice @cs, -$ex if $ex;
2012 8         43 @cs;
2013             }
2014              
2015             # Delete internal constants from namespace
2016             delete @PDL::DSP::Windows::{qw(
2017             I
2018             )};
2019              
2020             1;
2021              
2022             __END__