Processing math: 100%
LAL 7.7.0.1-6c6b863
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Modules Pages

Detailed Description

Performs complex-to-complex FFTs.

Synopsis

#include <lal/ComplexFFT.h>

Perform complex-to-complex fast Fourier transforms of vectors using the package FFTW [11] .

Description

This package provides a (X)LAL-style interface with the FFTW fast Fourier transform package [11] .

The routines LALCreateForwardComplexFFTPlan() and LALCreateReverseComplexFFTPlan() create plans for computing the forward and reverse FFTs of a given size. The optimum plan is either estimated (reasonably fast) if the measure flag is zero, or measured (can be time-consuming, but gives better performance) if the measure flag is non-zero. The routine LALDestroyComplexFFTPlan() destroys either of these flavours of plans.

The routine LALCOMPLEX8VectorFFT() performs either the forward or reverse FFT depending on the plan. The discrete Fourier transform H_k, k=0\ldots n-1 of a vector h_j, j=0\ldots n-1, of length n is defined by

H_k = \sum_{j=0}^{n-1} h_j e^{-2\pi ijk/n}

and, similarly, the inverse Fourier transform is defined by

h_j = \frac{1}{n}\sum_{k=0}^{n-1} H_k e^{2\pi ijk/n}.

However, the present implementation of the reverse FFT omits the factor of 1/n. The input and output vectors must be distinct.

Operating Instructions

const UINT4 n = 17;
ComplexFFTPlan *pfwd = NULL;
ComplexFFTPlan *prev = NULL;
COMPLEX8Vector *avec = NULL;
COMPLEX8Vector *bvec = NULL;
COMPLEX8Vector *cvec = NULL;
LALCreateForwardComplexFFTPlan( &status, &pfwd, n, 0 );
LALCreateReverseComplexFFTPlan( &status, &prev, n, 0 );
LALCCreateVector( &status, &avec, n );
LALCCreateVector( &status, &bvec, n );
LALCCreateVector( &status, &cvec, n );
<assign data>
LALCOMPLEX8VectorFFT( &status, bvec, avec, pfwd );
LALCOMPLEX8VectorFFT( &status, cvec, bvec, prev );
LALDestroyComplexFFTPlan( &status, &pfwd );
LALDestroyComplexFFTPlan( &status, &prev );
#define ComplexFFTPlan
Definition: ComplexFFT.h:83
uint32_t UINT4
Four-byte unsigned integer.
void LALCCreateVector(LALStatus *, COMPLEX8Vector **, UINT4)
void LALCDestroyVector(LALStatus *, COMPLEX8Vector **)
Vector of type COMPLEX8, see DATATYPE-Vector types for more details.
Definition: LALDatatypes.h:163
LAL status structure, see The LALStatus structure for more details.
Definition: LALDatatypes.h:947

Algorithm

The FFTW [11] is used.

Uses

Notes

  1. The sign convention used here is the opposite of the definition in Numerical Recipes [22] , but agrees with the one used by FFTW [11] and the other LIGO software components.
  2. The result of the inverse FFT must be multiplied by 1/n to recover the original vector. This is different from the datacondAPI where the factor is applied by default.
  3. The size n of the transform can be any positive integer; the performance is O(n\log n). However, better performance is obtained if n is the product of powers of 2, 3, 5, 7, and zero or one power of either 11 or 13. Transforms when n is a power of 2 are especially fast. See Ref. [11] .
  4. LALMalloc() is used by all the fftw routines.
  5. The input and output vectors for LALCOMPLEX8VectorFFT() must be distinct.

Prototypes

COMPLEX8FFTPlan * XLALCreateCOMPLEX8FFTPlan (UINT4 size, int fwdflg, int measurelvl)
 Returns a new COMPLEX8FFTPlan A COMPLEX8FFTPlan is required to perform a FFT that involves complex data. More...
 
COMPLEX8FFTPlan * XLALCreateForwardCOMPLEX8FFTPlan (UINT4 size, int measurelvl)
 Returns a new COMPLEX8FFTPlan for a forward transform. More...
 
COMPLEX8FFTPlan * XLALCreateReverseCOMPLEX8FFTPlan (UINT4 size, int measurelvl)
 Returns a new COMPLEX8FFTPlan for a reverse transform. More...
 
void XLALDestroyCOMPLEX8FFTPlan (COMPLEX8FFTPlan *plan)
 Destroys a COMPLEX8FFTPlan. More...
 
int XLALCOMPLEX8VectorFFT (COMPLEX8Vector *_LAL_RESTRICT_ output, const COMPLEX8Vector *_LAL_RESTRICT_ input, const COMPLEX8FFTPlan *plan)
 Perform a COMPLEX8Vector to COMPLEX8Vector FFT. More...
 
COMPLEX16FFTPlan * XLALCreateCOMPLEX16FFTPlan (UINT4 size, int fwdflg, int measurelvl)
 Returns a new COMPLEX16FFTPlan. More...
 
COMPLEX16FFTPlan * XLALCreateForwardCOMPLEX16FFTPlan (UINT4 size, int measurelvl)
 Returns a new COMPLEX16FFTPlan for a forward transform. More...
 
COMPLEX16FFTPlan * XLALCreateReverseCOMPLEX16FFTPlan (UINT4 size, int measurelvl)
 Returns a new COMPLEX16FFTPlan for a reverse transform. More...
 
void XLALDestroyCOMPLEX16FFTPlan (COMPLEX16FFTPlan *plan)
 Destroys a COMPLEX16FFTPlan. More...
 
int XLALCOMPLEX16VectorFFT (COMPLEX16Vector *_LAL_RESTRICT_ output, const COMPLEX16Vector *_LAL_RESTRICT_ input, const COMPLEX16FFTPlan *plan)
 Perform a COMPLEX16Vector to COMPLEX16Vector FFT. More...
 

Data Structures

struct  COMPLEX8FFTPlan
 Plan to perform an FFT of COMPLEX8 data. More...
 
struct  COMPLEX16FFTPlan
 Plan to perform an FFT of COMPLEX16 data. More...
 

Macros

#define SINGLE_PRECISION
 
#define tagComplexFFTPlan   tagCOMPLEX8FFTPlan
 
#define ComplexFFTPlan   COMPLEX8FFTPlan
 

Files

file  ComplexFFTTest.c
 Tests the routines in ComplexFFT.h.
 

Error Codes

#define COMPLEXFFTH_ENULL   1
 Null pointer. More...
 
#define COMPLEXFFTH_ENNUL   2
 Non-null pointer. More...
 
#define COMPLEXFFTH_ESIZE   4
 Invalid input size. More...
 
#define COMPLEXFFTH_ESZMM   8
 Size mismatch. More...
 
#define COMPLEXFFTH_ESLEN   16
 Invalid/mismatched sequence lengths. More...
 
#define COMPLEXFFTH_ESAME   32
 Input/Output data vectors are the same. More...
 
#define COMPLEXFFTH_EALOC   64
 Memory allocation failed. More...
 
#define COMPLEXFFTH_EFFTW   128
 Error in FFTW. More...
 
#define COMPLEXFFTH_ESNGL   256
 FFTW library is not single-precision. More...
 
#define COMPLEXFFTH_EINTL   512
 Error in Intel FFT library. More...
 
#define COMPLEXFFTH_ESIGN   1024
 Unknown sign of transform in plan. More...
 

Function Documentation

◆ XLALCreateCOMPLEX8FFTPlan()

COMPLEX8FFTPlan * XLALCreateCOMPLEX8FFTPlan ( UINT4  size,
int  fwdflg,
int  measurelvl 
)

Returns a new COMPLEX8FFTPlan A COMPLEX8FFTPlan is required to perform a FFT that involves complex data.

A different plan is required for each size of the complex data vectors and for each direction of transform (forward or reverse). A forward transform performs

Z[k] = \sum_{j=0}^{N-1} e^{-2\pi ijk/N}\,z[j]

where N, the size of the transform, is the length of the vectors z and Z. A reverse transform performs

w[j] = \sum_{k=0}^{N-1} e^{+2\pi ijk/N}\,W[k]

where N, the size of the transform, is the length of the vectors w and W.

Note
The reverse transform of the forward transform of some data is equal to N times the original data (we therefore call it a "reverse" transform rather than an "inverse" transform).
Parameters
[in]sizeThe number of points in the complex data.
[in]fwdflgSet non-zero for a forward FFT plan; otherwise create a reverse plan
[in]measurelvlMeasurement level for plan creation:
  • 0: no measurement, just estimate the plan;
  • 1: measure the best plan;
  • 2: perform a lengthy measurement of the best plan;
  • 3: perform an exhasutive measurement of the best plan.
Returns
A pointer to an allocated COMPLEX8FFTPlan structure is returned upon successful completion. Otherwise, a NULL pointer is returned and xlalErrno is set to indicate the error.
Errors:
The XLALCreateCOMPLEX8Plan() function shall fail if:
  • [XLAL_EBADLEN] The size of the requested plan is 0.
  • [XLAL_ENOMEM] Insufficient storage space is available.
  • [XLAL_EFAILED] The call to the underlying FFTW routine failed.

◆ XLALCreateForwardCOMPLEX8FFTPlan()

COMPLEX8FFTPlan * XLALCreateForwardCOMPLEX8FFTPlan ( UINT4  size,
int  measurelvl 
)

Returns a new COMPLEX8FFTPlan for a forward transform.

A COMPLEX8FFTPlan is required to perform a FFT that involves complex data. A different plan is required for each size of the complex data vectors. A forward transform performs

Z[k] = \sum_{j=0}^{N-1} e^{-2\pi ijk/N}\,z[j]

where N, the size of the transform, is the length of the vector z and Z.

Parameters
[in]sizeThe number of points in the complex data.
[in]measurelvlMeasurement level for plan creation:
  • 0: no measurement, just estimate the plan;
  • 1: measure the best plan;
  • 2: perform a lengthy measurement of the best plan;
  • 3: perform an exhasutive measurement of the best plan.
Returns
A pointer to an allocated COMPLEX8FFTPlan structure is returned upon successful completion. Otherwise, a NULL pointer is returned and xlalErrno is set to indicate the error.
Errors:
The XLALCreateForwardCOMPLEX8Plan() function shall fail if:
  • [XLAL_EBADLEN] The size of the requested plan is 0.
  • [XLAL_ENOMEM] Insufficient storage space is available.
  • [XLAL_EFAILED] The call to the underlying FFTW routine failed.

Definition at line 94 of file CudaComplexFFT.c.

◆ XLALCreateReverseCOMPLEX8FFTPlan()

COMPLEX8FFTPlan * XLALCreateReverseCOMPLEX8FFTPlan ( UINT4  size,
int  measurelvl 
)

Returns a new COMPLEX8FFTPlan for a reverse transform.

A COMPLEX8FFTPlan is required to perform a FFT that involves complex data. A reverse transform performs

w[j] = \sum_{k=0}^{N-1} e^{+2\pi ijk/N}\,W[k]

where N, the size of the transform, is the length of the vectors w and W.

Note
The reverse transform of the forward transform of some data is equal to N times the original data (we therefore call it a "reverse" transform rather than an "inverse" transform).
Parameters
[in]sizeThe number of points in the complex data.
[in]measurelvlMeasurement level for plan creation:
  • 0: no measurement, just estimate the plan;
  • 1: measure the best plan;
  • 2: perform a lengthy measurement of the best plan;
  • 3: perform an exhasutive measurement of the best plan.
Returns
A pointer to an allocated COMPLEX8FFTPlan structure is returned upon successful completion. Otherwise, a NULL pointer is returned and xlalErrno is set to indicate the error.
Errors:
The XLALCreateReverseCOMPLEX8Plan() function shall fail if:
  • [XLAL_EBADLEN] The size of the requested plan is 0.
  • [XLAL_ENOMEM] Insufficient storage space is available.
  • [XLAL_EFAILED] The call to the underlying FFTW routine failed.

Definition at line 104 of file CudaComplexFFT.c.

◆ XLALDestroyCOMPLEX8FFTPlan()

void XLALDestroyCOMPLEX8FFTPlan ( COMPLEX8FFTPlan *  plan)

Destroys a COMPLEX8FFTPlan.

Parameters
[in]planA pointer to the COMPLEX8FFTPlan to be destroyed.

Definition at line 114 of file CudaComplexFFT.c.

◆ XLALCOMPLEX8VectorFFT()

int XLALCOMPLEX8VectorFFT ( COMPLEX8Vector *_LAL_RESTRICT_  output,
const COMPLEX8Vector *_LAL_RESTRICT_  input,
const COMPLEX8FFTPlan *  plan 
)

Perform a COMPLEX8Vector to COMPLEX8Vector FFT.

This routine computes

Z[k] = \sum_{j=0}^{N-1} e^{\mp2\pi ijk/N}\,z[j],

and where the minus sign is used if a forward plan is provided as the argument and the plus sign is used if a reverse plan is provided as the argument; here N is the length of the input and output vectors z and Z.

Parameters
[out]outputThe complex output data vector Z of length N
[in]inputThe input complex data vector z of length N
[in]planThe FFT plan to use for the transform
Note
The input and output vectors must be distinct.
Returns
0 upon successful completion or non-zero upon failure.
Errors:
The XLALCOMPLEX8VectorFFT() function shall fail if:
  • [XLAL_EFAULT] A NULL pointer is provided as one of the arguments.
  • [XLAL_EINVAL] A argument is invalid or the input and output data vectors are the same.
  • [XLAL_EBADLEN] The input vector, output vector, and plan size are incompatible.
  • [XLAL_ENOMEM] Insufficient storage space is available.

Definition at line 133 of file CudaComplexFFT.c.

◆ XLALCreateCOMPLEX16FFTPlan()

COMPLEX16FFTPlan * XLALCreateCOMPLEX16FFTPlan ( UINT4  size,
int  fwdflg,
int  measurelvl 
)

Returns a new COMPLEX16FFTPlan.

A COMPLEX16FFTPlan is required to perform a FFT that involves complex data. A different plan is required for each size of the complex data vectors and for each direction of transform (forward or reverse). A forward transform performs

Z[k] = \sum_{j=0}^{N-1} e^{-2\pi ijk/N}\,z[j]

where N, the size of the transform, is the length of the vectors z and Z. A reverse transform performs

w[j] = \sum_{k=0}^{N-1} e^{+2\pi ijk/N}\,W[k]

where N, the size of the transform, is the length of the vectors w and W.

Note
The reverse transform of the forward transform of some data is equal to N times the original data (we therefore call it a "reverse" transform rather than an "inverse" transform).
Parameters
[in]sizeThe number of points in the complex data.
[in]fwdflgSet non-zero for a forward FFT plan; otherwise create a reverse plan
[in]measurelvlMeasurement level for plan creation:
  • 0: no measurement, just estimate the plan;
  • 1: measure the best plan;
  • 2: perform a lengthy measurement of the best plan;
  • 3: perform an exhasutive measurement of the best plan.
Returns
A pointer to an allocated COMPLEX16FFTPlan structure is returned upon successful completion. Otherwise, a NULL pointer is returned and xlalErrno is set to indicate the error.
Errors:
The XLALCreateCOMPLEX16Plan() function shall fail if:
  • [XLAL_EBADLEN] The size of the requested plan is 0.
  • [XLAL_ENOMEM] Insufficient storage space is available.
  • [XLAL_EFAILED] The call to the underlying FFTW routine failed.

Definition at line 169 of file CudaComplexFFT.c.

◆ XLALCreateForwardCOMPLEX16FFTPlan()

COMPLEX16FFTPlan * XLALCreateForwardCOMPLEX16FFTPlan ( UINT4  size,
int  measurelvl 
)

Returns a new COMPLEX16FFTPlan for a forward transform.

A COMPLEX16FFTPlan is required to perform a FFT that involves complex data. A different plan is required for each size of the complex data vectors. A forward transform performs

Z[k] = \sum_{j=0}^{N-1} e^{-2\pi ijk/N}\,z[j]

where N, the size of the transform, is the length of the vector z and Z.

Parameters
[in]sizeThe number of points in the complex data.
[in]measurelvlMeasurement level for plan creation:
  • 0: no measurement, just estimate the plan;
  • 1: measure the best plan;
  • 2: perform a lengthy measurement of the best plan;
  • 3: perform an exhasutive measurement of the best plan.
Returns
A pointer to an allocated COMPLEX16FFTPlan structure is returned upon successful completion. Otherwise, a NULL pointer is returned and xlalErrno is set to indicate the error.
Errors:
The XLALCreateForwardCOMPLEX16Plan() function shall fail if:
  • [XLAL_EBADLEN] The size of the requested plan is 0.
  • [XLAL_ENOMEM] Insufficient storage space is available.
  • [XLAL_EFAILED] The call to the underlying FFTW routine failed.

Definition at line 236 of file CudaComplexFFT.c.

◆ XLALCreateReverseCOMPLEX16FFTPlan()

COMPLEX16FFTPlan * XLALCreateReverseCOMPLEX16FFTPlan ( UINT4  size,
int  measurelvl 
)

Returns a new COMPLEX16FFTPlan for a reverse transform.

A COMPLEX16FFTPlan is required to perform a FFT that involves complex data. A reverse transform performs

w[j] = \sum_{k=0}^{N-1} e^{+2\pi ijk/N}\,W[k]

where N, the size of the transform, is the length of the vectors w and W.

Note
The reverse transform of the forward transform of some data is equal to N times the original data (we therefore call it a "reverse" transform rather than an "inverse" transform).
Parameters
[in]sizeThe number of points in the complex data.
[in]measurelvlMeasurement level for plan creation:
  • 0: no measurement, just estimate the plan;
  • 1: measure the best plan;
  • 2: perform a lengthy measurement of the best plan;
  • 3: perform an exhasutive measurement of the best plan.
Returns
A pointer to an allocated COMPLEX16FFTPlan structure is returned upon successful completion. Otherwise, a NULL pointer is returned and xlalErrno is set to indicate the error.
Errors:
The XLALCreateReverseCOMPLEX16Plan() function shall fail if:
  • [XLAL_EBADLEN] The size of the requested plan is 0.
  • [XLAL_ENOMEM] Insufficient storage space is available.
  • [XLAL_EFAILED] The call to the underlying FFTW routine failed.

Definition at line 246 of file CudaComplexFFT.c.

◆ XLALDestroyCOMPLEX16FFTPlan()

void XLALDestroyCOMPLEX16FFTPlan ( COMPLEX16FFTPlan *  plan)

Destroys a COMPLEX16FFTPlan.

Parameters
[in]planA pointer to the COMPLEX16FFTPlan to be destroyed.

Definition at line 256 of file CudaComplexFFT.c.

◆ XLALCOMPLEX16VectorFFT()

int XLALCOMPLEX16VectorFFT ( COMPLEX16Vector *_LAL_RESTRICT_  output,
const COMPLEX16Vector *_LAL_RESTRICT_  input,
const COMPLEX16FFTPlan *  plan 
)

Perform a COMPLEX16Vector to COMPLEX16Vector FFT.

This routine computes

Z[k] = \sum_{j=0}^{N-1} e^{\mp2\pi ijk/N}\,z[j],

and where the minus sign is used if a forward plan is provided as the argument and the plus sign is used if a reverse plan is provided as the argument; here N is the length of the input and output vectors z and Z.

Parameters
[out]outputThe complex output data vector Z of length N
[in]inputThe input complex data vector z of length N
[in]planThe FFT plan to use for the transform
Note
The input and output vectors must be distinct.
Returns
0 upon successful completion or non-zero upon failure.
Errors:
The XLALCOMPLEX16VectorFFT() function shall fail if:
  • [XLAL_EFAULT] A NULL pointer is provided as one of the arguments.
  • [XLAL_EINVAL] A argument is invalid or the input and output data vectors are the same.
  • [XLAL_EBADLEN] The input vector, output vector, and plan size are incompatible.
  • [XLAL_ENOMEM] Insufficient storage space is available.

Definition at line 273 of file CudaComplexFFT.c.

Macro Definition Documentation

◆ SINGLE_PRECISION

#define SINGLE_PRECISION

Definition at line 145 of file ComplexFFT.c.

◆ COMPLEXFFTH_ENULL

#define COMPLEXFFTH_ENULL   1

Null pointer.

Definition at line 51 of file ComplexFFT.h.

◆ COMPLEXFFTH_ENNUL

#define COMPLEXFFTH_ENNUL   2

Non-null pointer.

Definition at line 52 of file ComplexFFT.h.

◆ COMPLEXFFTH_ESIZE

#define COMPLEXFFTH_ESIZE   4

Invalid input size.

Definition at line 53 of file ComplexFFT.h.

◆ COMPLEXFFTH_ESZMM

#define COMPLEXFFTH_ESZMM   8

Size mismatch.

Definition at line 54 of file ComplexFFT.h.

◆ COMPLEXFFTH_ESLEN

#define COMPLEXFFTH_ESLEN   16

Invalid/mismatched sequence lengths.

Definition at line 55 of file ComplexFFT.h.

◆ COMPLEXFFTH_ESAME

#define COMPLEXFFTH_ESAME   32

Input/Output data vectors are the same.

Definition at line 56 of file ComplexFFT.h.

◆ COMPLEXFFTH_EALOC

#define COMPLEXFFTH_EALOC   64

Memory allocation failed.

Definition at line 57 of file ComplexFFT.h.

◆ COMPLEXFFTH_EFFTW

#define COMPLEXFFTH_EFFTW   128

Error in FFTW.

Definition at line 58 of file ComplexFFT.h.

◆ COMPLEXFFTH_ESNGL

#define COMPLEXFFTH_ESNGL   256

FFTW library is not single-precision.

Definition at line 59 of file ComplexFFT.h.

◆ COMPLEXFFTH_EINTL

#define COMPLEXFFTH_EINTL   512

Error in Intel FFT library.

Definition at line 60 of file ComplexFFT.h.

◆ COMPLEXFFTH_ESIGN

#define COMPLEXFFTH_ESIGN   1024

Unknown sign of transform in plan.

Definition at line 61 of file ComplexFFT.h.

◆ tagComplexFFTPlan

#define tagComplexFFTPlan   tagCOMPLEX8FFTPlan

Definition at line 82 of file ComplexFFT.h.

◆ ComplexFFTPlan

#define ComplexFFTPlan   COMPLEX8FFTPlan

Definition at line 83 of file ComplexFFT.h.