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
CudaRealFFT.c
Go to the documentation of this file.
1/*
2* Copyright (C) 2010 Karsten Wiesner, Jolien Creighton, Kipp Cannon, Shin Kee Chung
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#include <string.h>
21
22#include <lal/LALDatatypes.h>
23#include <lal/LALMalloc.h>
24#include <lal/RealFFT.h>
25#include <lal/XLALError.h>
26#include <lal/FFTWMutex.h>
27#include <fftw3.h>
28
29#include "CudaFunctions.h"
30#include "CudaFFT.h"
31
32#ifdef __GNUC__
33#define UNUSED __attribute__ ((unused))
34#else
35#define UNUSED
36#endif
37
38struct
40{
43 cufftHandle plan;
46};
47
48struct
50{
53 fftw_plan plan;
54};
55
56
57/*
58 *
59 * REAL4 XLAL Functions
60 *
61 */
62
63
64REAL4FFTPlan * XLALCreateREAL4FFTPlan( UINT4 size, int fwdflg, UNUSED int measurelvl )
65{
66 UINT4 createSize;
67 REAL4FFTPlan *plan;
68
69 if ( ! size )
71
72 /* "use" measurelvl */
73 measurelvl = 0;
74
75 /* allocate memory for the plan and the temporary arrays */
76 plan = XLALMalloc( sizeof( *plan ) );
77 if ( ! plan )
78 {
79 XLALFree( plan );
81 }
82
83 /*
84 * Use a different size for plan creation to avoid current CUDA bug
85 * in performing FFT to array with size 1
86 */
87
88 if( size == 1 ) createSize = 2;
89 else createSize = size;
90 /* LAL_FFTW_WISDOM_LOCK; */
91 if ( fwdflg ) /* forward */
92 cufftPlan1d( &plan->plan, createSize, CUFFT_R2C, 1 );
93 else /* reverse */
94 cufftPlan1d( &plan->plan, createSize, CUFFT_C2R, 1 );
95 /* LAL_FFTW_WISDOM_UNLOCK; */
96 /* check to see success of plan creation */
97
98 /* "Plan=0" Bugfix by Wiesner, K.: plan->plan is an integer handle not a pointer and 0 is a valid handle
99 So checking against 0 and occasionaly destroy the plan is a bug.
100
101 if ( ! plan->plan )
102 {
103 XLALFree( plan );
104 XLAL_ERROR_NULL( XLAL_EFAILED );
105 }
106 */
107
108 /* Allocate memory in the GPU */
109 plan->d_real = XLALCudaMallocReal(size);
110 plan->d_complex = XLALCudaMallocComplex(size/2 + 1);
111
112 /* now set remaining plan fields */
113 plan->size = size;
114 plan->sign = ( fwdflg ? -1 : 1 );
115
116 return plan;
117}
118
119
120REAL4FFTPlan * XLALCreateForwardREAL4FFTPlan( UINT4 size, int measurelvl )
121{
122 REAL4FFTPlan *plan;
123 plan = XLALCreateREAL4FFTPlan( size, 1, measurelvl );
124 if ( ! plan )
126 return plan;
127}
128
129
130REAL4FFTPlan * XLALCreateReverseREAL4FFTPlan( UINT4 size, int measurelvl )
131{
132 REAL4FFTPlan *plan;
133 plan = XLALCreateREAL4FFTPlan( size, 0, measurelvl );
134 if ( ! plan )
136 return plan;
137}
138
139
140void XLALDestroyREAL4FFTPlan( REAL4FFTPlan *plan )
141{
142 if ( ! plan )
144
145 /* Plan=0 Bugfix
146 if ( ! plan->plan )
147 XLAL_ERROR_VOID( XLAL_EINVAL );
148 */
149
150 /* LAL_FFTW_WISDOM_LOCK; */
151 /* Free the Cuda specific variables */
152 XLALCudaFree( plan->d_real );
153 cufftDestroy( plan->plan );
154 XLALCudaFree( plan->d_complex );
155 /* LAL_FFTW_WISDOM_UNLOCK; */
156 memset( plan, 0, sizeof( *plan ) );
157 XLALFree( plan );
158 return;
159}
160
161int XLALREAL4ForwardFFT( COMPLEX8Vector *output, const REAL4Vector *input, const REAL4FFTPlan *plan )
162{
163 if ( ! output || ! input || ! plan )
165 /* Plan=0 Bugfix
166 if ( ! plan->plan || ! plan->size || plan->sign != -1 )
167 */
168 if ( ! plan->size || plan->sign != -1 )
170 if ( ! output->data || ! input->data )
172 if ( input->length != plan->size || output->length != plan->size/2 + 1 )
174
175 /* do the fft */
176 /*
177 * Special check for vector with size 1 to avoid CUDA bug
178 */
179 if( plan->size == 1 )
180 {
181 output->data[0] = crectf( input->data[0], 0.0 );
182 }
183 else
184 cudafft_execute_r2c( plan->plan,
185 (cufftComplex *)(output->data), (cufftReal *)(input->data),
186 (cufftComplex *)plan->d_complex, (cufftReal *)plan->d_real, plan->size );
187
188 /* Nyquist frequency */
189 if( plan->size%2 == 0 )
190 {
191 output->data[plan->size/2] = crectf( crealf(output->data[plan->size/2]), 0.0 );
192 }
193
194 return 0;
195}
196
197
198int XLALREAL4ReverseFFT( REAL4Vector *output, const COMPLEX8Vector *input, const REAL4FFTPlan *plan )
199{
200 if ( ! output || ! input || ! plan )
202 /* Plan=0 Bugfix
203 if ( ! plan->plan || ! plan->size || plan->sign != 1 )
204 */
205 if ( ! plan->size || plan->sign != 1 )
207 if ( ! output->data || ! input->data )
209 if ( output->length != plan->size || input->length != plan->size/2 + 1 )
211 if ( cimagf(input->data[0]) != 0.0 )
212 XLAL_ERROR( XLAL_EDOM ); /* imaginary part of DC must be zero */
213 if ( ! plan->size % 2 && cimagf(input->data[plan->size/2]) != 0.0 )
214 XLAL_ERROR( XLAL_EDOM ); /* imaginary part of Nyquist must be zero */
215
216 /* perform the fft */
217 /*
218 * Special check to handle array of size 1
219 */
220 if( plan->size == 1 )
221 {
222 output->data[0] = crealf(input->data[0]);
223 }
224 else
225 cudafft_execute_c2r( plan->plan,
226 (cufftReal *)(output->data), (cufftComplex *)(input->data),
227 (cufftReal *)plan->d_real, (cufftComplex *)plan->d_complex, plan->size );
228
229 return 0;
230}
231
232
233int XLALREAL4VectorFFT( REAL4Vector * _LAL_RESTRICT_ output, const REAL4Vector * _LAL_RESTRICT_ input, const REAL4FFTPlan *plan )
234{
235 COMPLEX8 *tmp;
236 UINT4 k;
237
238 if ( ! output || ! input || ! plan )
240 /* Plan=0 Bugfix
241 if ( ! plan->plan || ! plan->size )
242 */
243 if ( ! plan->size )
245 if ( ! output->data || ! input->data || output->data == input->data )
246 XLAL_ERROR( XLAL_EINVAL ); /* note: must be out-of-place */
247 if ( output->length != plan->size || input->length != plan->size )
249
250 if( plan->size == 1 )
251 {
252 output->data[0] = input->data[0];
253 }
254 else
255 {
256 tmp = XLALMalloc( plan->size * sizeof(*tmp) );
257 /* do the fft */
258 if( plan->sign == -1 )
259 {
260 cudafft_execute_r2c( plan->plan,
261 (cufftComplex *)tmp, (cufftReal *)(input->data),
262 (cufftComplex *)plan->d_complex, (cufftReal *)plan->d_real, plan->size );
263 output->data[0] = crealf(tmp[0]);
264
265 for( k = 1; k < (plan->size + 1)/2; k++ )
266 {
267 output->data[k] = crealf(tmp[k]);
268 output->data[plan->size - k] = cimagf(tmp[k]);
269 }
270
271 if( plan->size % 2 == 0 )
272 {
273 output->data[plan->size/2] = crealf(tmp[plan->size/2]);
274 }
275 }
276 else
277 {
278 tmp[0] = crectf( input->data[0], 0.0 );
279
280 for( k = 1; k < (plan->size + 1)/2; k++ )
281 {
282 tmp[k] = crectf( input->data[k], input->data[plan->size - k] );
283 }
284
285 if( plan->size%2 == 0 )
286 {
287 tmp[plan->size/2] = crectf( input->data[plan->size/2], 0.0 );
288 }
289
290 cudafft_execute_c2r( plan->plan,
291 (cufftReal *)(output->data), (cufftComplex *)tmp,
292 (cufftReal *)plan->d_real, (cufftComplex *)plan->d_complex, plan->size );
293 }
294 XLALFree(tmp);
295 }
296
297 return 0;
298}
299
300
302 const REAL4FFTPlan *plan )
303{
304 COMPLEX8 *tmp;
305 UINT4 k;
306
307 if ( ! spec || ! data || ! plan )
309 /* Plan=0 Bugfix
310 if ( ! plan->plan || ! plan->size )
311 */
312 if (! plan->size )
314 if ( ! spec->data || ! data->data )
316 if ( data->length != plan->size || spec->length != plan->size/2 + 1 )
318
319 /* allocate temporary storage space */
320 tmp = XLALMalloc( (plan->size/2 + 1) * sizeof( *tmp ) );
321 if ( ! tmp )
323
324 /* Check for size 1 to avoid the CUDA bug */
325 if( plan->size == 1 )
326 tmp[0] = crectf( data->data[0], cimagf(tmp[0]) );
327 /* transform the data */
328 else
329 cudafft_execute_r2c( plan->plan,
330 (cufftComplex *)(tmp), (cufftReal *)(data->data),
331 (cufftComplex *)plan->d_complex, (cufftReal *)plan->d_real, plan->size );
332
333 /* now reconstruct the spectrum from the temporary storage */
334
335 /* dc component */
336 spec->data[0] = crealf(tmp[0]) * crealf(tmp[0]);
337
338 /* other components */
339 for (k = 1; k < (plan->size + 1)/2; ++k) /* k < size/2 rounded up */
340 {
341 REAL4 re = crealf(tmp[k]);
342 REAL4 im = cimagf(tmp[k]);
343 spec->data[k] = re * re + im * im;
344 spec->data[k] *= 2.0; /* accounts for negative frequency part */
345 }
346
347 /* Nyquist frequency */
348 if ( plan->size%2 == 0 ) /* size is even */
349 spec->data[k] = crealf(tmp[k]) * crealf(tmp[k]);
350
351 /* clenup and exit */
352 XLALFree( tmp );
353 return 0;
354}
355
356
357
358/*
359 *
360 * REAL8 XLAL Functions
361 *
362 */
363
364
365
366REAL8FFTPlan * XLALCreateREAL8FFTPlan( UINT4 size, int fwdflg, int measurelvl )
367{
368 REAL8FFTPlan *plan;
369 REAL8 *tmp1;
370 REAL8 *tmp2;
371 int flags = FFTW_UNALIGNED;
372
373 if ( ! size )
375
376 /* based on measurement level, set fftw3 flags to perform
377 * requested degree of measurement */
378 switch ( measurelvl )
379 {
380 case 0: /* estimate */
381 flags |= FFTW_ESTIMATE;
382 break;
383 default: /* exhaustive measurement */
384 flags |= FFTW_EXHAUSTIVE;
385 /* fall-through */
386 case 2: /* lengthy measurement */
387 flags |= FFTW_PATIENT;
388 /* fall-through */
389 case 1: /* measure the best plan */
390 flags |= FFTW_MEASURE;
391 break;
392 }
393
394 /* allocate memory for the plan and the temporary arrays */
395 plan = XLALMalloc( sizeof( *plan ) );
396 tmp1 = XLALMalloc( size * sizeof( *tmp1 ) );
397 tmp2 = XLALMalloc( size * sizeof( *tmp2 ) );
398 if ( ! plan || ! tmp1 || ! tmp2 )
399 {
400 XLALFree( plan );
401 XLALFree( tmp1 );
402 XLALFree( tmp2 );
404 }
405
407 if ( fwdflg ) /* forward */
408 plan->plan = fftw_plan_r2r_1d( size, tmp1, tmp2, FFTW_R2HC, flags );
409 else /* reverse */
410 plan->plan = fftw_plan_r2r_1d( size, tmp1, tmp2, FFTW_HC2R, flags );
412
413 /* free temporary arrays */
414 XLALFree( tmp2 );
415 XLALFree( tmp1 );
416
417 /* check to see success of plan creation */
418 /* Plan=0 Bugfix
419 if ( ! plan->plan )
420 {
421 XLALFree( plan );
422 XLAL_ERROR_NULL( XLAL_EFAILED );
423 }
424 */
425
426 /* now set remaining plan fields */
427 plan->size = size;
428 plan->sign = ( fwdflg ? -1 : 1 );
429
430 return plan;
431}
432
433
434REAL8FFTPlan * XLALCreateForwardREAL8FFTPlan( UINT4 size, int measurelvl )
435{
436 REAL8FFTPlan *plan;
437 plan = XLALCreateREAL8FFTPlan( size, 1, measurelvl );
438 if ( ! plan )
440 return plan;
441}
442
443
444REAL8FFTPlan * XLALCreateReverseREAL8FFTPlan( UINT4 size, int measurelvl )
445{
446 REAL8FFTPlan *plan;
447 plan = XLALCreateREAL8FFTPlan( size, 0, measurelvl );
448 if ( ! plan )
450 return plan;
451}
452
453
454void XLALDestroyREAL8FFTPlan( REAL8FFTPlan *plan )
455{
456 if ( ! plan )
458 /* Plan=0 Bugfix
459 if ( ! plan->plan )
460 XLAL_ERROR_VOID( XLAL_EINVAL );
461 */
463 fftw_destroy_plan( plan->plan );
465 memset( plan, 0, sizeof( *plan ) );
466 XLALFree( plan );
467
468 return;
469}
470
471int XLALREAL8ForwardFFT( COMPLEX16Vector *output, const REAL8Vector *input, const REAL8FFTPlan *plan )
472{
473 REAL8 *tmp;
474 UINT4 k;
475
476 if ( ! output || ! input || ! plan )
478 /* Plan=0 Bugfix
479 if ( ! plan->plan || ! plan->size || plan->sign != -1 )
480 */
481 if ( ! plan->size || plan->sign != -1 )
483 if ( ! output->data || ! input->data )
485 if ( input->length != plan->size || output->length != plan->size/2 + 1 )
487
488 /* create temporary storage space */
489 tmp = XLALMalloc( plan->size * sizeof( *tmp ) );
490 if ( ! tmp )
492
493 /* do the fft */
494 fftw_execute_r2r( plan->plan, input->data, tmp );
495
496 /* now unpack the results into the output vector */
497
498 /* dc component */
499 output->data[0] = crect( tmp[0], 0.0 );
500
501 /* other components */
502 for ( k = 1; k < (plan->size + 1)/2; ++k ) /* k < size/2 rounded up */
503 {
504 output->data[k] = crect( tmp[k], tmp[plan->size - k] );
505 }
506
507 /* Nyquist frequency */
508 if ( plan->size%2 == 0 ) /* n is even */
509 {
510 output->data[plan->size/2] = crect( tmp[plan->size/2], 0.0 );
511 }
512
513 XLALFree( tmp );
514 return 0;
515}
516
517
518int XLALREAL8ReverseFFT( REAL8Vector *output, const COMPLEX16Vector *input, const REAL8FFTPlan *plan )
519{
520 REAL8 *tmp;
521 UINT4 k;
522
523 if ( ! output || ! input || ! plan )
525 /* Plan=0 Bugfix:
526 if ( ! plan->plan || ! plan->size || plan->sign != 1 )
527 */
528 if ( ! plan->size || plan->sign != 1 )
530 if ( ! output->data || ! input->data )
532 if ( output->length != plan->size || input->length != plan->size/2 + 1 )
534 if ( cimag(input->data[0]) != 0.0 )
535 XLAL_ERROR( XLAL_EDOM ); /* imaginary part of DC must be zero */
536 if ( ! plan->size % 2 && cimag(input->data[plan->size/2]) != 0.0 )
537 XLAL_ERROR( XLAL_EDOM ); /* imaginary part of Nyquist must be zero */
538
539 /* create temporary storage space */
540 tmp = XLALMalloc( plan->size * sizeof( *tmp ) );
541 if ( ! tmp )
543
544 /* unpack input into temporary array */
545
546 /* dc component */
547 tmp[0] = creal(input->data[0]);
548
549 /* other components */
550 for ( k = 1; k < (plan->size + 1)/2; ++k ) /* k < size / 2 rounded up */
551 {
552 tmp[k] = creal(input->data[k]);
553 tmp[plan->size - k] = cimag(input->data[k]);
554 }
555
556 /* Nyquist component */
557 if ( plan->size%2 == 0 ) /* n is even */
558 tmp[plan->size/2] = creal(input->data[plan->size/2]);
559
560 /* perform the fft */
561 fftw_execute_r2r( plan->plan, tmp, output->data );
562
563 /* cleanup and exit */
564 XLALFree( tmp );
565 return 0;
566}
567
568
569int XLALREAL8VectorFFT( REAL8Vector * _LAL_RESTRICT_ output, const REAL8Vector * _LAL_RESTRICT_ input, const REAL8FFTPlan *plan )
570{
571 if ( ! output || ! input || ! plan )
573 /* Plan=0 Bugfix
574 if ( ! plan->plan || ! plan->size )
575 */
576 if (! plan->size )
578 if ( ! output->data || ! input->data || output->data == input->data )
579 XLAL_ERROR( XLAL_EINVAL ); /* note: must be out-of-place */
580 if ( output->length != plan->size || input->length != plan->size )
582
583 /* do the fft */
584 fftw_execute_r2r( plan->plan, input->data, output->data );
585 return 0;
586}
587
588
590 const REAL8FFTPlan *plan )
591{
592 REAL8 *tmp;
593 UINT4 k;
594
595 if ( ! spec || ! data || ! plan )
597 /* Plan=0 Bugfix
598 if ( ! plan->plan || ! plan->size )
599 */
600 if ( ! plan->size )
602 if ( ! spec->data || ! data->data )
604 if ( data->length != plan->size || spec->length != plan->size/2 + 1 )
606
607 /* allocate temporary storage space */
608 tmp = XLALMalloc( plan->size * sizeof( *tmp ) );
609 if ( ! tmp )
611
612 /* transform the data */
613 fftw_execute_r2r( plan->plan, data->data, tmp );
614
615 /* now reconstruct the spectrum from the temporary storage */
616
617 /* dc component */
618 spec->data[0] = tmp[0] * tmp[0];
619
620 /* other components */
621 for (k = 1; k < (plan->size + 1)/2; ++k) /* k < size/2 rounded up */
622 {
623 REAL8 re = tmp[k];
624 REAL8 im = tmp[plan->size - k];
625 spec->data[k] = re * re + im * im;
626 spec->data[k] *= 2.0; /* accounts for negative frequency part */
627 }
628
629 /* Nyquist frequency */
630 if ( plan->size%2 == 0 ) /* size is even */
631 spec->data[plan->size/2] = tmp[plan->size/2] * tmp[plan->size/2];
632
633 /* clenup and exit */
634 XLALFree( tmp );
635 return 0;
636}
int cudafft_execute_c2r(cufftHandle plan, cufftReal *output, const cufftComplex *input, cufftReal *d_output, cufftComplex *d_input, UINT4 size)
Definition: CudaFFT.cu:66
int cudafft_execute_r2c(cufftHandle plan, cufftComplex *output, const cufftReal *input, cufftComplex *d_output, cufftReal *d_input, UINT4 size)
Definition: CudaFFT.cu:50
REAL4 * XLALCudaMallocReal(UINT4 size)
Definition: CudaFunctions.c:7
void XLALCudaFree(void *d_data)
Definition: CudaFunctions.c:29
COMPLEX8 * XLALCudaMallocComplex(UINT4 size)
Definition: CudaFunctions.c:18
int XLALREAL4PowerSpectrum(REAL4Vector *spec, const REAL4Vector *data, const REAL4FFTPlan *plan)
Definition: CudaRealFFT.c:301
int XLALREAL8VectorFFT(REAL8Vector *_LAL_RESTRICT_ output, const REAL8Vector *_LAL_RESTRICT_ input, const REAL8FFTPlan *plan)
Definition: CudaRealFFT.c:569
REAL4FFTPlan * XLALCreateREAL4FFTPlan(UINT4 size, int fwdflg, UNUSED int measurelvl)
Definition: CudaRealFFT.c:64
#define LAL_FFTW_WISDOM_UNLOCK
Definition: FFTWMutex.h:37
#define LAL_FFTW_WISDOM_LOCK
Definition: FFTWMutex.h:36
#define _LAL_RESTRICT_
Definition: LALStddef.h:35
#define crect(re, im)
Construct a COMPLEX16 from real and imaginary parts.
double REAL8
Double precision real floating-point number (8 bytes).
#define crectf(re, im)
Construct a COMPLEX8 from real and imaginary parts.
uint32_t UINT4
Four-byte unsigned integer.
float complex COMPLEX8
Single-precision floating-point complex number (8 bytes total)
int32_t INT4
Four-byte signed integer.
float REAL4
Single precision real floating-point number (4 bytes).
#define XLALMalloc(n)
Definition: LALMalloc.h:44
#define XLALFree(p)
Definition: LALMalloc.h:47
int XLALREAL8ForwardFFT(COMPLEX16Vector *output, const REAL8Vector *input, const REAL8FFTPlan *plan)
Performs a forward FFT of REAL8 data.
Definition: CudaRealFFT.c:471
int XLALREAL8ReverseFFT(REAL8Vector *output, const COMPLEX16Vector *input, const REAL8FFTPlan *plan)
Performs a reverse FFT of REAL8 data.
Definition: CudaRealFFT.c:518
void XLALDestroyREAL8FFTPlan(REAL8FFTPlan *plan)
Destroys a REAL8FFTPlan.
Definition: CudaRealFFT.c:454
REAL4FFTPlan * XLALCreateReverseREAL4FFTPlan(UINT4 size, int measurelvl)
Returns a new REAL4FFTPlan for a reverse transform.
Definition: CudaRealFFT.c:130
REAL4FFTPlan * XLALCreateForwardREAL4FFTPlan(UINT4 size, int measurelvl)
Returns a new REAL4FFTPlan for a forward transform.
Definition: CudaRealFFT.c:120
int XLALREAL8PowerSpectrum(REAL8Vector *spec, const REAL8Vector *data, const REAL8FFTPlan *plan)
Computes the power spectrum of REAL8 data.
Definition: CudaRealFFT.c:589
void XLALDestroyREAL4FFTPlan(REAL4FFTPlan *plan)
Destroys a REAL4FFTPlan.
Definition: CudaRealFFT.c:140
REAL8FFTPlan * XLALCreateReverseREAL8FFTPlan(UINT4 size, int measurelvl)
Returns a new REAL8FFTPlan for a reverse transform.
Definition: CudaRealFFT.c:444
int XLALREAL4VectorFFT(REAL4Vector *_LAL_RESTRICT_ output, const REAL4Vector *_LAL_RESTRICT_ input, const REAL4FFTPlan *plan)
Perform a REAL4Vector to REAL4Vector FFT.
Definition: CudaRealFFT.c:233
REAL8FFTPlan * XLALCreateForwardREAL8FFTPlan(UINT4 size, int measurelvl)
Returns a new REAL8FFTPlan for a forward transform.
Definition: CudaRealFFT.c:434
int XLALREAL4ReverseFFT(REAL4Vector *output, const COMPLEX8Vector *input, const REAL4FFTPlan *plan)
Performs a reverse FFT of REAL4 data.
Definition: CudaRealFFT.c:198
REAL8FFTPlan * XLALCreateREAL8FFTPlan(UINT4 size, int fwdflg, int measurelvl)
Returns a new REAL8FFTPlan.
Definition: CudaRealFFT.c:366
int XLALREAL4ForwardFFT(COMPLEX8Vector *output, const REAL4Vector *input, const REAL4FFTPlan *plan)
Performs a forward FFT of REAL4 data.
Definition: CudaRealFFT.c:161
#define XLAL_ERROR_VOID(...)
Macro to invoke a failure from a XLAL routine returning void.
Definition: XLALError.h:726
#define XLAL_ERROR_NULL(...)
Macro to invoke a failure from a XLAL routine returning a pointer.
Definition: XLALError.h:713
#define XLAL_ERROR(...)
Macro to invoke a failure from a XLAL routine returning an integer.
Definition: XLALError.h:700
@ XLAL_EBADLEN
Inconsistent or invalid length.
Definition: XLALError.h:419
@ XLAL_ENOMEM
Memory allocation error.
Definition: XLALError.h:407
@ XLAL_EFAULT
Invalid pointer.
Definition: XLALError.h:408
@ XLAL_EFUNC
Internal function call failed bit: "or" this with existing error number.
Definition: XLALError.h:462
@ XLAL_EDOM
Input domain error.
Definition: XLALError.h:410
@ XLAL_EINVAL
Invalid argument.
Definition: XLALError.h:409
Vector of type COMPLEX16, see DATATYPE-Vector types for more details.
Definition: LALDatatypes.h:172
UINT4 length
Number of elements in array.
Definition: LALDatatypes.h:176
COMPLEX16 * data
Pointer to the data array.
Definition: LALDatatypes.h:177
Vector of type COMPLEX8, see DATATYPE-Vector types for more details.
Definition: LALDatatypes.h:163
UINT4 length
Number of elements in array.
Definition: LALDatatypes.h:167
COMPLEX8 * data
Pointer to the data array.
Definition: LALDatatypes.h:168
Vector of type REAL4, see DATATYPE-Vector types for more details.
Definition: LALDatatypes.h:145
REAL4 * data
Pointer to the data array.
Definition: LALDatatypes.h:150
UINT4 length
Number of elements in array.
Definition: LALDatatypes.h:149
Vector of type REAL8, see DATATYPE-Vector types for more details.
Definition: LALDatatypes.h:154
REAL8 * data
Pointer to the data array.
Definition: LALDatatypes.h:159
UINT4 length
Number of elements in array.
Definition: LALDatatypes.h:158
Plan to perform FFT of REAL4 data.
Definition: CudaRealFFT.c:40
cufftHandle plan
Definition: CudaRealFFT.c:43
INT4 sign
sign in transform exponential, -1 for forward, +1 for reverse
Definition: CudaRealFFT.c:41
COMPLEX8 * d_complex
Definition: CudaRealFFT.c:45
UINT4 size
length of the real data vector for this plan
Definition: CudaRealFFT.c:42
REAL4 * d_real
Definition: CudaRealFFT.c:44
Plan to perform FFT of REAL8 data.
Definition: CudaRealFFT.c:50
UINT4 size
length of the real data vector for this plan
Definition: CudaRealFFT.c:52
fftw_plan plan
the FFTW plan
Definition: CudaRealFFT.c:53
INT4 sign
sign in transform exponential, -1 for forward, +1 for reverse
Definition: CudaRealFFT.c:51
void output(int gps_sec, int output_type)
Definition: tconvert.c:440