Loading [MathJax]/extensions/TeX/AMSsymbols.js
LAL 7.7.0.1-5e288d3
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
Window.h
Go to the documentation of this file.
1/*
2 * Copyright (C) 2007 Bruce Allen, Duncan Brown, Jolien Creighton, Kipp
3 * Cannon, Teviet Creighton
4 *
5 * This program is free software; you can redistribute it and/or modify it
6 * under the terms of the GNU General Public License as published by the
7 * Free Software Foundation; either version 2 of the License, or (at your
8 * option) any later version.
9 *
10 * This program is distributed in the hope that it will be useful, but
11 * WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13 * General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License along
16 * with with program; see the file COPYING. If not, write to the Free
17 * Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
18 * 02110-1301 USA
19 */
20
21#ifndef _WINDOW_H
22#define _WINDOW_H
23
24#include <lal/LALDatatypes.h>
25
26#ifdef __cplusplus
27extern "C" {
28#endif
29
30/**
31 * \defgroup Window_h Header Window.h
32 * \ingroup lal_window
33 * \brief This header file provides routines and structures to create and store window functions (also called a taper,
34 * lag window, or apodization function).
35 *
36 * ### Synopsis ###
37 *
38 * \code
39 * #include <lal/Window.h>
40 * \endcode
41 *
42 * The \ref Window_h package defines two data types, the \c REAL4Window and
43 * \c REAL8Window types. These are suitable for storing window function
44 * data, providing storage for a sequence of samples as well as metadata
45 * about the window such as the sum-of-squarse of the samples. The package
46 * also provides a collection of functions to construct a variety of
47 * pre-defined window functions.
48 *
49 * The use of window functions in signal analysis is well documented in
50 * many places.
51 *
52 * Definitions of most of the window functions can be found in Numerical
53 * Recipes \cite ptvf1992 equations 13.4.13-13.4.15. Definitions of the
54 * remaining windows can be found in <em>Spectral analysis for physical
55 * applications</em> \cite pw93 Section 6.11. Definition of the Kaiser
56 * window can be found in <em>Discrete-time Signal Processing</em> by
57 * Oppenheim and Schafer, p.474.
58 * *
59 * ### Description ###
60 *
61 * These functions create or destroy a time-domain window function in a vector
62 * of specified length. If you wish to construct a custom window, call
63 * <tt>XLALCreateRectangularREAL8Window()</tt> (or the \c REAL4
64 * version), then replace the samples inside it with your own, and update the
65 * \c sumofsquares and \c sum elements. If the window function
66 * proves useful, consider adding it here so that others can benefit.
67 *
68 * It is convenient to describe the windows as functions on the normalized
69 * domain \f$y \in [-1, 1]\f$. The window is zero outside this domain. The
70 * window functions defined in this package are as follows.
71 *
72 * ### Rectangle ###
73 *
74 * \f{equation}{
75 * w(y)
76 * = 1.
77 * \f}</dd>
78 *
79 * ### Hann ###
80 *
81 * \f{equation}{
82 * w(y)
83 * = \cos^2 \frac{\pi}{2} y.
84 * \f}</dd>
85 *
86 * ### Welch ###
87 *
88 * \f{equation}{
89 * w(y)
90 * = 1 - y^2.
91 * \f}</dd>
92 *
93 * ### Bartlett ###
94 *
95 * \f{equation}{
96 * w(y)
97 * = 1 - |y|.
98 * \f}</dd>
99 *
100 * ### Parzen ###
101 *
102 * \f{equation}{
103 * w(y)
104 * = \left\{ \begin{array}{ll}
105 * 1 - 6 y^2 (1 - |y|) & |y| \leq 1 / 2, \\
106 * 2 (1 - |y|)^3 & |y| > 1 / 2.
107 * \end{array}\right.
108 * \f}</dd>
109 *
110 * ### Papoulis ###
111 *
112 * \f{equation}{
113 * w(y)
114 * = \frac{1}{\pi} \sin \pi |y| + (1 - |y|) \cos \pi |y|.
115 * \f}</dd>
116 *
117 * ### Hamming ###
118 *
119 * \f{equation}{
120 * w(y)
121 * = 0.08 + 0.92 \cos^{2} \frac{\pi}{2} y.
122 * \f}
123 * This is the same as the Hann window, but with an additional DC bias, or
124 * &quot;foot,&quot; of 0.08.</dd>
125 *
126 * ### Kaiser ###
127 *
128 * \f{equation}{
129 * w(y)
130 * = I_0 \left( \beta \sqrt{1-y^2} \right) / I_0(\beta),
131 * \f}
132 * where \f$I_0(x)\f$ is the \f$0\f$th order, modified Bessel function of the first
133 * kind. The shape parameter \f$\beta \in [0, \infty]\f$ sets the sharpness of
134 * the central peak. \f$\beta = 0\f$ yields the rectangle window, \f$\beta
135 * \rightarrow \infty\f$ yields a \f$\delta\f$ function with a single non-zero
136 * sample in the middle. This window is difficult to compute for large
137 * \f$\beta\f$, and an asymptotic approximation is used for \f$\beta \ge 705\f$. A
138 * linearly-interpolated transition occurs between \f$\beta = 695\f$ and \f$\beta =
139 * 705\f$. Finite-difference derivatives of the window with respect to \f$\beta\f$
140 * are unlikely to work well in this regime.</dd>
141 *
142 * ### Creighton ###
143 *
144 * \f{equation}{
145 * w(y)
146 * = \exp \left( -\beta \frac{y^2}{1 - y^2} \right).
147 * \f}
148 * This window function is based on a fairly standard \f$C_{\infty}\f$ test
149 * function used in distribution theory, e.g.\ <em>Green's Functions and
150 * Boundary Value Problems</em> \cite stakgold79 , by Stakgold. The shape parameter
151 * \f$\beta \in [0, \infty]\f$ sets the sharpness of the central peak. \f$\beta =
152 * 0\f$ yields the rectangle window, \f$\beta \rightarrow \infty\f$ yields a
153 * \f$\delta\f$ function with a single non-zero sample in the middle.</dd>
154 *
155 * ### Tukey ###
156 *
157 * \f{equation}{
158 * w(y)
159 * = \left\{ \begin{array}{ll}
160 * \sin^2 \left[ \frac{\pi}{2} (|y| - 1) / \beta \right] & |y| \geq 1 -
161 * \beta, \\
162 * 1 & |y| < 1 - \beta.
163 * \end{array} \right.
164 * \f}
165 * The shape parameter \f$\beta \in [0, 1]\f$ sets what fraction of the window is
166 * spanned by the tapers. \f$\beta = 0\f$ yields the rectangle window, \f$\beta =
167 * 1\f$ yields the Hann window.</dd>
168 *
169 * ### Gauss ###
170 *
171 * \f{equation}{
172 * w(y)
173 * = \exp \left( -\frac{1}{2} \beta^{2} y^{2} \right).
174 * \f}
175 * The shape parameter \f$\beta \in [0, \infty]\f$ sets the sharpness of the
176 * central peak. \f$\beta = 0\f$ yields the rectangle window, \f$\beta \rightarrow \infty\f$ yields a \f$\delta\f$ function with a single non-zero sample in the
177 * middle.
178 *
179 * ### Lanczos ###
180 *
181 * \f{equation}{
182 * w(y)
183 * = \frac{\sin \pi y}{\pi y}.
184 * \f}
185 * The Lanczos window is the central lobe of the sinc function. This is
186 * used, for example, in finite impulse response resampling to modulate a
187 * truncated sinc interpolating kernel.
188 *
189 * These window functions are shown in \ref window_t "this figure", showing various windows as functions of the normalized
190 * independend variable \f$y\f$, choosing \f$\beta = 6\f$ for the Kaiser window, \f$\beta = 2\f$ for the Creighton window,
191 * \f$\beta = 0.5\f$ for the Tukey window, and \f$\beta = 3\f$ for the Gauss window.
192 *
193 * \anchor window_t
194 * \image html window_t.png "Various windows as functions of the normalized independend variable y"
195 *
196 * For a vector of length \f$L\f$ (an integer), the mapping from integer array
197 * index \f$i\f$ to normalized co-ordinate \f$y\f$ is
198 * \f{equation}{
199 * y(i)
200 * = \left\{ \begin{array}{ll}
201 * 0 & L \le 1, \\
202 * 2 i / (L - 1) - 1 & L > 1,
203 * \end{array} \right.
204 * \f}
205 * where \f$0 \le i < L\f$, and floating-point division is used. This agrees with
206 * J. G. Proakis and D. G. Manolakis, <em>Digital Signal Processing</em>
207 * \cite pm95 , and \c MatLab. The first sample is \f$y = -1\f$, the last
208 * sample is \f$y = +1\f$. For odd-lengthed vectors, the middle sample is \f$y =
209 * 0\f$, while for even-lengthed vectors \f$y = 0\f$ occurs half-way between the two
210 * middle samples. Substituting \f$y(i)\f$ into the definitions of the window
211 * functions above yields \f$w(i)\f$, the value of the window function at the
212 * integer sample \f$i\f$.
213 *
214 * The Fourier transforms of the windows are shown as functions of \f$1 / y\f$ in
215 * \ref window_f "this figure", showing frequency behaviour of various windows as functions
216 * of the inverse of the normalized independend variable \f$y\f$, choosing \f$\beta = 6\f$
217 * for the Kaiser window, \f$\beta = 2\f$ for the Creighton window, \f$\beta = 0.5\f$ for
218 * the Tukey window, and \f$\beta = 3\f$ for the Gauss window.
219 *
220 * \anchor window_f
221 * \image html window_f.png "Frequency behaviour of various windows as functions of the inverse of the normalized independend variable y"
222 *
223 * Since the Fourier transform of windowed data is the Fourier transform of
224 * the data convolved with the Fourier transform of the window,
225 * \ref window_f "this figure" is the major guideline for selecting a window. One
226 * can see that windows with a narrow central lobe tend to have higher
227 * sidelobes, and windows which suppress their low-order sidelobes tend to
228 * have more power in the high-order sidelobes. The choice of window thus
229 * depends on whether one is trying to resolve nearby spectral features of
230 * comparable magnitude (suggesting a rectangular or a Welch window), to
231 * reduce spectral bias and low-order sidelobes (a Hamming or Kaiser window),
232 * or to measure a broad spectrum with a large dynamical range (a Creighton or
233 * a Papoulis window).
234 *
235 */
236/** @{ */
237
238
239/**
240 * Structure for storing REAL4 window function data, providing storage for a sequence of samples
241 * as well as metadata about the window such as the sum-of-squarse of the samples
242 */
243typedef struct tagREAL4Window {
244 REAL4Sequence *data; /**< The window function samples */
245 REAL8 sumofsquares; /**< The sum of the squares of the window function samples */
246 REAL8 sum; /**< The sum of the window function samples */
248
249
250/**
251 * Structure for storing REAL8 window function data, providing storage for a sequence of samples
252 * as well as metadata about the window such as the sum-of-squarse of the samples
253 */
254typedef struct tagREAL8Window {
255 REAL8Sequence *data; /**< The window function samples */
256 REAL8 sumofsquares; /**< The sum of the squares of the window function samples */
257 REAL8 sum; /**< The sum of the window function samples */
259
260
261#ifdef SWIG /* SWIG interface directives */
262SWIGLAL(OWNS_THIS_ARG(REAL4Sequence*, sequence));
263SWIGLAL(OWNS_THIS_ARG(REAL8Sequence*, sequence));
264#endif
267#ifdef SWIG /* SWIG interface directives */
268SWIGLAL_CLEAR(OWNS_THIS_ARG(REAL4Sequence*, sequence));
269SWIGLAL_CLEAR(OWNS_THIS_ARG(REAL8Sequence*, sequence));
270#endif
271
284
297
300
305
306int XLALCheckNamedWindow ( const char *windowName, const BOOLEAN haveBeta );
307
308REAL8Window *XLALCreateNamedREAL8Window ( const char *windowName, REAL8 beta, UINT4 length );
309REAL4Window *XLALCreateNamedREAL4Window ( const char *windowName, REAL8 beta, UINT4 length );
310
311/** @} */
312
313#ifdef __cplusplus
314}
315#endif
316
317#endif /* _WINDOW_H */
unsigned char BOOLEAN
Boolean logical type, see Headers LAL(Atomic)Datatypes.h for more details.
double REAL8
Double precision real floating-point number (8 bytes).
uint32_t UINT4
Four-byte unsigned integer.
float REAL4
Single precision real floating-point number (4 bytes).
COMPLEX8Sequence * XLALUnitaryWindowCOMPLEX8Sequence(COMPLEX8Sequence *sequence, const REAL4Window *window)
Single-precision complex version of XLALUnitaryWindowREAL8Sequence().
Definition: Window.c:324
REAL4Window * XLALCreateParzenREAL4Window(UINT4 length)
Definition: Window.c:709
REAL8Window * XLALCreateBartlettREAL8Window(UINT4 length)
Definition: Window.c:402
REAL8Window * XLALCreatePapoulisREAL8Window(UINT4 length)
Definition: Window.c:442
REAL8Window * XLALCreateGaussREAL8Window(UINT4 length, REAL8 beta)
Definition: Window.c:623
REAL4Window * XLALCreateREAL4WindowFromSequence(REAL4Sequence *sequence)
Single-precision version of XLALCreateREAL8WindowFromSequence().
Definition: Window.c:225
REAL8Window * XLALCreateNamedREAL8Window(const char *windowName, REAL8 beta, UINT4 length)
Generic window-function wrapper, allowing to select a window by its name.
Definition: Window.c:877
REAL4Window * XLALCreateHammingREAL4Window(UINT4 length)
Definition: Window.c:721
REAL4Window * XLALCreateLanczosREAL4Window(UINT4 length)
Definition: Window.c:751
REAL8Window * XLALCreateLanczosREAL8Window(UINT4 length)
Definition: Window.c:650
REAL4Window * XLALCreatePapoulisREAL4Window(UINT4 length)
Definition: Window.c:715
REAL8Sequence * XLALUnitaryWindowREAL8Sequence(REAL8Sequence *sequence, const REAL8Window *window)
Multiply a REAL8Sequence in-place by a REAL8Window with a normalization that preserves the variance o...
Definition: Window.c:264
REAL8Window * XLALCreateREAL8WindowFromSequence(REAL8Sequence *sequence)
Constructs a new REAL8Window from a REAL8Sequence.
Definition: Window.c:204
REAL8Window * XLALCreateCreightonREAL8Window(UINT4 length, REAL8 beta)
Definition: Window.c:563
COMPLEX16Sequence * XLALUnitaryWindowCOMPLEX16Sequence(COMPLEX16Sequence *sequence, const REAL8Window *window)
Double-precision complex version of XLALUnitaryWindowREAL8Sequence().
Definition: Window.c:284
REAL8Window * XLALCreateHammingREAL8Window(UINT4 length)
Definition: Window.c:461
REAL8Window * XLALCreateWelchREAL8Window(UINT4 length)
Definition: Window.c:384
REAL4Window * XLALCreateTukeyREAL4Window(UINT4 length, REAL4 beta)
Definition: Window.c:739
REAL8Window * XLALCreateTukeyREAL8Window(UINT4 length, REAL8 beta)
Definition: Window.c:597
void XLALDestroyREAL8Window(REAL8Window *window)
Definition: Window.c:668
REAL4Window * XLALCreateKaiserREAL4Window(UINT4 length, REAL4 beta)
Definition: Window.c:727
REAL8Window * XLALCreateHannREAL8Window(UINT4 length)
Definition: Window.c:367
REAL4Window * XLALCreateNamedREAL4Window(const char *windowName, REAL8 beta, UINT4 length)
Definition: Window.c:936
int XLALCheckNamedWindow(const char *windowName, const BOOLEAN haveBeta)
Check whether a named window-function exists and whether it requires a parameter.
Definition: Window.c:854
void XLALDestroyREAL4Window(REAL4Window *window)
Definition: Window.c:757
REAL8Window * XLALCreateParzenREAL8Window(UINT4 length)
Definition: Window.c:421
REAL8Window * XLALCreateRectangularREAL8Window(UINT4 length)
Definition: Window.c:350
REAL4Window * XLALCreateCreightonREAL4Window(UINT4 length, REAL4 beta)
Definition: Window.c:733
REAL4Window * XLALCreateGaussREAL4Window(UINT4 length, REAL4 beta)
Definition: Window.c:745
REAL4Window * XLALCreateHannREAL4Window(UINT4 length)
Definition: Window.c:691
REAL4Sequence * XLALUnitaryWindowREAL4Sequence(REAL4Sequence *sequence, const REAL4Window *window)
Single-precision version of XLALUnitaryWindowREAL8Sequence().
Definition: Window.c:304
REAL4Window * XLALCreateRectangularREAL4Window(UINT4 length)
Definition: Window.c:685
REAL4Window * XLALCreateBartlettREAL4Window(UINT4 length)
Definition: Window.c:703
REAL4Window * XLALCreateWelchREAL4Window(UINT4 length)
Definition: Window.c:697
REAL8Window * XLALCreateKaiserREAL8Window(UINT4 length, REAL8 beta)
Definition: Window.c:478
Vector of type COMPLEX16, see DATATYPE-Vector types for more details.
Definition: LALDatatypes.h:172
Vector of type COMPLEX8, see DATATYPE-Vector types for more details.
Definition: LALDatatypes.h:163
Vector of type REAL4, see DATATYPE-Vector types for more details.
Definition: LALDatatypes.h:145
Structure for storing REAL4 window function data, providing storage for a sequence of samples as well...
Definition: Window.h:243
REAL4Sequence * data
The window function samples.
Definition: Window.h:244
REAL8 sum
The sum of the window function samples.
Definition: Window.h:246
REAL8 sumofsquares
The sum of the squares of the window function samples.
Definition: Window.h:245
Vector of type REAL8, see DATATYPE-Vector types for more details.
Definition: LALDatatypes.h:154
Structure for storing REAL8 window function data, providing storage for a sequence of samples as well...
Definition: Window.h:254
REAL8Sequence * data
The window function samples.
Definition: Window.h:255
REAL8 sumofsquares
The sum of the squares of the window function samples.
Definition: Window.h:256
REAL8 sum
The sum of the window function samples.
Definition: Window.h:257