Loading [MathJax]/extensions/TeX/AMSsymbols.js
LALPulsar 7.1.1.1-b246709
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages
genSFTWindows.c
Go to the documentation of this file.
1/*
2 * Copyright (C) 2022 Karl Wette
3 *
4 * This program is free software; you can redistribute it and/or modify
5 * it under the terms of the GNU General Public License as published by
6 * the Free Software Foundation; either version 2 of the License, or
7 * (at your option) any later version.
8 *
9 * This program is distributed in the hope that it will be useful,
10 * but WITHOUT ANY WARRANTY; without even the implied warranty of
11 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 * GNU General Public License for more details.
13 *
14 * You should have received a copy of the GNU General Public License
15 * along with with program; see the file COPYING. If not, write to the
16 * Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston,
17 * MA 02110-1301 USA
18 */
19
20/*
21 * Author: K. Wette
22 *
23 * A simple script to print out SFT window functions
24 */
25
26#include "config.h"
27
28#include <math.h>
29
30#include <lal/LALStdlib.h>
31#include <lal/AVFactories.h>
32#include <lal/Window.h>
33
34/* GLOBAL VARIABLES */
35extern REAL8 winFncRMS;
38
39/* FUNCTION PROTOTYPES */
40int WindowData( REAL8 r );
41int WindowDataTukey2( void );
42int WindowDataHann( void );
43
44int main( void )
45{
46
47 const size_t NWINDOWS = 4;
48 const size_t WINDOWLENGTH = 256 * 1800;
49
50 // default command line arguments from MakeSFTs.c
51 const REAL8 windowR = 0.001;
52
53 // allocate memory
54 char windownames[NWINDOWS][1024];
55 REAL8Vector *windows[NWINDOWS];
57 for ( size_t i = 0; i < NWINDOWS; ++i ) {
58 windows[i] = XLALCreateREAL8Vector( WINDOWLENGTH );
59 XLAL_CHECK_MAIN( windows[i] != NULL, XLAL_ENOMEM );
60 for ( size_t j = 0; j < WINDOWLENGTH; ++j ) {
61 windows[i]->data[j] = 1.0;
62 }
63 }
64
65 size_t w = 0;
66
67 {
68 snprintf( windownames[w], sizeof( windownames[w] ), "lalpulsar_MakeSFTs Matlab style Tukey window [windowR=%g]", windowR );
69 dataDouble.data = windows[w];
70 WindowData( windowR );
72 }
73
74 ++w;
75
76 {
77 snprintf( windownames[w], sizeof( windownames[w] ), "lalpulsar_MakeSFTs Hann window" );
78 dataDouble.data = windows[w];
81 }
82
83 ++w;
84
85 {
86 REAL8 param = windowR;
87 REAL8Window *win = XLALCreateTukeyREAL8Window( WINDOWLENGTH, param );
88 XLAL_CHECK_MAIN( win != NULL, XLAL_EFUNC );
89 snprintf( windownames[w], sizeof( windownames[w] ), "XLALCreateTukeyREAL8Window(param=%g)", param );
90 for ( size_t j = 0; j < WINDOWLENGTH; ++j ) {
91 windows[w]->data[j] *= win->data->data[j];
92 }
93 winrms[w] = sqrt( win->sumofsquares / win->data->length );
95 }
96
97 ++w;
98
99 {
100 REAL8Window *win = XLALCreateHannREAL8Window( WINDOWLENGTH );
101 XLAL_CHECK_MAIN( win != NULL, XLAL_EFUNC );
102 snprintf( windownames[w], sizeof( windownames[w] ), "XLALCreateHannREAL8Window()" );
103 for ( size_t j = 0; j < WINDOWLENGTH; ++j ) {
104 windows[w]->data[j] *= win->data->data[j];
105 }
106 winrms[w] = sqrt( win->sumofsquares / win->data->length );
108 }
109
111
112 // output windows
113 for ( size_t i = 0; i < NWINDOWS; ++i ) {
114 printf( "%s%c", windownames[i], i + 1 < NWINDOWS ? ',' : '\n' );
115 }
116 for ( size_t i = 0; i < NWINDOWS; ++i ) {
117 printf( "%0.8f%c", winrms[i], i + 1 < NWINDOWS ? ',' : '\n' );
118 }
119 for ( size_t j = 0; j < WINDOWLENGTH; ++j ) {
120 for ( size_t i = 0; i < NWINDOWS; ++i ) {
121 printf( "%0.8f%c", windows[i]->data[j], i + 1 < NWINDOWS ? ',' : '\n' );
122 }
123 }
124
125 // cleanup
126 for ( size_t i = 0; i < NWINDOWS; ++i ) {
127 XLALDestroyREAL8Vector( windows[i] );
128 }
130
131 return XLAL_SUCCESS;
132
133}
int j
void LALCheckMemoryLeaks(void)
#define NWINDOWS
const double w
REAL4TimeSeries dataSingle
REAL8 winFncRMS
int WindowDataTukey2(void)
int WindowData(REAL8 r)
REAL8TimeSeries dataDouble
int main(void)
Definition: genSFTWindows.c:44
int WindowDataHann(void)
double REAL8
static const INT4 r
REAL8Vector * XLALCreateREAL8Vector(UINT4 length)
void XLALDestroyREAL8Vector(REAL8Vector *vector)
REAL8Window * XLALCreateTukeyREAL8Window(UINT4 length, REAL8 beta)
void XLALDestroyREAL8Window(REAL8Window *window)
REAL8Window * XLALCreateHannREAL8Window(UINT4 length)
#define XLAL_CHECK_MAIN(assertion,...)
XLAL_ENOMEM
XLAL_SUCCESS
XLAL_EFUNC
XLAL_EFAILED
float data[BLOCKSIZE]
Definition: hwinject.c:360
REAL8Sequence * data
REAL8 * data