LAL  7.5.0.1-08ee4f4

Detailed Description

Performs real-to-complex and complex-to-real FFTs.

Synopsis

#include <lal/RealFFT.h>

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

XLAL Functions

Synopsis

#include <lal/RealFFT.h>
REAL4FFTPlan * XLALCreateREAL4FFTPlan( UINT4 size, int fwdflg, int measurelvl );
REAL4FFTPlan * XLALCreateForwardREAL4FFTPlan( UINT4 size, int measurelvl );
REAL4FFTPlan * XLALCreateReverseREAL4FFTPlan( UINT4 size, int measurelvl );
void XLALDestroyREAL4FFTPlan( REAL4FFTPlan *plan );
int XLALREAL4ForwardFFT( COMPLEX8Vector *output, REAL4Vector *input, REAL4FFTPlan *plan );
int XLALREAL4ReverseFFT( REAL4Vector *output, COMPLEX8Vector *input, REAL4FFTPlan *plan );
int XLALREAL4VectorFFT( REAL4Vector *output, REAL4Vector *input, REAL4FFTPlan *plan );
int XLALREAL4PowerSpectrum( REAL4Vector *spec, REAL4Vector *data, REAL4FFTPlan *plan );
REAL8FFTPlan * XLALCreateREAL8FFTPlan( UINT4 size, int fwdflg, int measurelvl );
REAL8FFTPlan * XLALCreateForwardREAL8FFTPlan( UINT4 size, int measurelvl );
REAL8FFTPlan * XLALCreateReverseREAL8FFTPlan( UINT4 size, int measurelvl );
void XLALDestroyREAL8FFTPlan( REAL8FFTPlan *plan );
int XLALREAL8ForwardFFT( COMPLEX16Vector *output, REAL8Vector *input, REAL8FFTPlan *plan );
int XLALREAL8ReverseFFT( REAL8Vector *output, COMPLEX16Vector *input, REAL8FFTPlan *plan );
int XLALREAL8VectorFFT( REAL8Vector *output, REAL8Vector *input, REAL8FFTPlan *plan );
int XLALREAL8PowerSpectrum( REAL8Vector *spec, REAL8Vector *data, REAL8FFTPlan *plan );
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
uint32_t UINT4
Four-byte unsigned integer.
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
REAL4FFTPlan * XLALCreateReverseREAL4FFTPlan(UINT4 size, int measurelvl)
Returns a new REAL4FFTPlan for a reverse transform.
Definition: CudaRealFFT.c:130
void XLALDestroyREAL8FFTPlan(REAL8FFTPlan *plan)
Destroys a REAL8FFTPlan.
Definition: CudaRealFFT.c:454
REAL4FFTPlan * XLALCreateForwardREAL4FFTPlan(UINT4 size, int measurelvl)
Returns a new REAL4FFTPlan for a forward transform.
Definition: CudaRealFFT.c:120
REAL8FFTPlan * XLALCreateReverseREAL8FFTPlan(UINT4 size, int measurelvl)
Returns a new REAL8FFTPlan for a reverse transform.
Definition: CudaRealFFT.c:444
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 * XLALCreateREAL8FFTPlan(UINT4 size, int fwdflg, int measurelvl)
Returns a new REAL8FFTPlan.
Definition: CudaRealFFT.c:366
int XLALREAL4VectorFFT(REAL4Vector *_LAL_RESTRICT_ output, const REAL4Vector *_LAL_RESTRICT_ input, const REAL4FFTPlan *plan)
Perform a REAL4Vector to REAL4Vector FFT.
Definition: CudaRealFFT.c:233
int XLALREAL4ReverseFFT(REAL4Vector *output, const COMPLEX8Vector *input, const REAL4FFTPlan *plan)
Performs a reverse FFT of REAL4 data.
Definition: CudaRealFFT.c:198
REAL8FFTPlan * XLALCreateForwardREAL8FFTPlan(UINT4 size, int measurelvl)
Returns a new REAL8FFTPlan for a forward transform.
Definition: CudaRealFFT.c:434
int XLALREAL4ForwardFFT(COMPLEX8Vector *output, const REAL4Vector *input, const REAL4FFTPlan *plan)
Performs a forward FFT of REAL4 data.
Definition: CudaRealFFT.c:161
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
Vector of type REAL8, see DATATYPE-Vector types for more details.
Definition: LALDatatypes.h:154
void output(int gps_sec, int output_type)
Definition: tconvert.c:440

Description

The REAL4 routines are described below. These use single-precision FFTs, i.e., they convert REAL4Vectors into COMPLEX8Vectors and vice-versa. The REAL8 versions of the routines are the same but they are double-precision versions, i.e., they convert REAL8Vectors into COMPLEX16Vectors.

The routine XLALCreateREAL4FFTPlan() creates a REAL4FFTPlan structure to perform FFTs of vectors of length size. If fwdflg is non-zero then the plan is created to perform forward (real-to-complex) FFTs with a negative exponential sign. Otherwise the plan is created to perform reverse (complex-to-real) FFTs with a positive exponential sign. The value of measurelvl determines how much optimization of the plan FFTW will do with the most optimization taking the most amount of time. Reasonable values for measurelvl would be 0 for the fasted plan creation (FFTW does not measure the speed of any transform with this level but rather estimates which plan will be the fastet) or 1 to measure a few likely plans to determine the fastest.

XLALCreateForwardREAL4FFTPlan() is equivalent to XLALCreateREAL4FFTPlan() with fwdflg set to 1. XLALCreateReverseREAL4FFTPlan() is equivalent to XLALCreateREAL4FFTPlan() with fwdflg set to 0.

XLALDestroyREAL4FFTPlan() is used to destroy the plan, freeing all memory that was allocated in the structure as well as the structure itself. It can be used on either forward or reverse plans.

XLALREAL4ForwardFFT() and XLALREAL4ReverseFFT() perform forward (real to complex) and reverse (complex to real) transforms respectively. The plan supplied to these routines must be correctly generated for the direction of the transform. I.e., XLALREAL4ForwardFFT() cannot be supplied with a plan generated by XLALCreateReverseREAL4FFTPlan().

XLALREAL4VectorFFT() is a low-level routine that transforms a real vector to a half-complex real vector (with a forward plan) or a half-complex real vector to a real vector (with a reverse plan). If you're not sure what this means, don't use this routine. The input and output vectors (and their data) must be distinct pointers.

XLALREAL4PowerSpectrum() computes a real power spectrum of the input real vector and a forward FFT plan.

Return Values

Upon success, XLALCreateREAL4FFTPlan(), XLALCreateForwardREAL4FFTPlan(), and XLALCreateReverseREAL4FFTPlan() return a pointer to a newly-allocated FFT plan. Upon failure, they return a NULL pointer and set xlalErrno to one of the following values: XLAL_EBADLEN if size is not greater than zero, XLAL_ENOMEM if a memory allocation failed, or XLAL_EFAILED if the FFTW plan creation routine failed.

XLALDestroyREAL4FFTPlan() does not return any value but, upon failure, it will set xlalErrno to one of the following values: XLAL_EFAULT if the routine is provided a NULL pointer, or XLAL_EINVAL if the contents of the plan are invalid (e.g., if the routine is provided a plan that had been previously destroyed).

XLALREAL4ForwardFFT(), XLALREAL4ReverseFFT(), XLALREAL4VectorFFT(), and XLALREAL4PowerSpectrum() return the value 0 upon succes; upon failure they return XLAL_FAILURE and set xlalErrno to one of the following values: XLAL_EFAULT if one of the input pointers is NULL, XLAL_EINVAL if the input, output, or plan structures appears invalid or if the routine is passed a plan for the wrong transform directions or if the input and output data pointers are not distinct for XLALREAL4VectorFFT(), XLAL_EBADLEN if the input and output vectors and the plan have incompatible lengths, XLAL_ENOMEM if a memory allocation of temporary internal memory fails.

As before, the REAL8 versions of these routines behave the same way but for double-precision transforms.

LAL-style functions [DEPRECATED]

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

The routines LALCreateForwardRealFFTPlan() and LALCreateReverseRealFFTPlan() create plans for computing the forward (real-to-complex) and reverse (complex-to-real) FFTs of a specified 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 LALDestroyRealFFTPlan() destroys any of these flavours of plans.

The routines LALForwardRealFFT() and LALReverseRealFFT() perform the forward (real-to-complex) and reverse (complex-to-real) FFTs using the plans. The discrete Fourier transform \(H_k\), \(k=0\ldots\lfloor{n/2}\rfloor\) ( \(n/2\) rounded down), 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} \]

where \(H_k\) for \(\lfloor{n/2}\rfloor<k<n\) can be obtained from the relation \(H_k=H_{n-k}^\ast\). The present implementation of the reverse FFT omits the factor of \(1/n\).

The routines in this package require that the vector \(h_j\), \(j=0\ldots n-1\) be real; consequently, \(H_k=H_{n-k}^\ast\) ( \(0\le k\le\lfloor n/2\rfloor\)), i.e., the negative frequency Fourier components are the complex conjugate of the positive frequency Fourier components when the data is real. Therefore, one need compute and store only the first \(\lfloor n/2\rfloor+1\) components of \(H_k\); only the values of \(H_k\) for \(k=0\ldots \lfloor n/2\rfloor\) are returned (integer division is rounded down, e.g., \(\lfloor 7/2\rfloor=3\)).

The routine LALRealPowerSpectrum() computes the power spectrum \(P_k=2|H_k|^2\), \(k=1\ldots \lfloor (n-1)/2\rfloor\), \(P_0=|H_0|^2\), and \(P_{n/2}=|H_{n/2}|^2\) if \(n\) is even, of the data \(h_j\), \(j=0\ldots n-1\). The factor of two except at DC and Nyquist accounts for the power in negative frequencies.

The routine LALREAL4VectorFFT() is essentially a direct calls to FFTW routines without any re-packing of the data. This routine should not be used unless the user understands the packing used in FFTW.

Operating Instructions

const UINT4 n = 32;
RealFFTPlan *pfwd = NULL;
RealFFTPlan *prev = NULL;
REAL4Vector *hvec = NULL;
COMPLEX8Vector *Hvec = NULL;
REAL4Vector *Pvec = NULL;
LALCreateForwardRealFFTPlan( &status, &pfwd, n );
LALCreateReverseRealFFTPlan( &status, &prev, n );
LALSCreateVector( &status, &hvec, n );
LALCCreateVector( &status, &Hvec, n/2 + 1 );
LALSCreateVector( &status, &Pvec, n/2 + 1 );
<assign data>
LALRealPowerSpectrum( &status, Pvec, hvec, pfwd );
LALForwardRealFFT( &status, Hvec, hvec, pfwd );
LALReverseRealFFT( &status, hvec, Hvec, pinv );
LALDestroyRealFFTPlan( &status, &pfwd );
LALDestroyRealFFTPlan( &status, &pinv );
#define RealFFTPlan
Definition: RealFFT.h:195
void LALCCreateVector(LALStatus *, COMPLEX8Vector **, UINT4)
void LALCDestroyVector(LALStatus *, COMPLEX8Vector **)
void LALSDestroyVector(LALStatus *, REAL4Vector **)
void LALSCreateVector(LALStatus *, REAL4Vector **, UINT4)
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 Numerical Recipes [22] , but agrees with the one used by FFTW [11] and the other LIGO software components.
  2. The result of the reverse FFT must be multiplied by \(1/n\) to recover the original vector. This is unlike the Numerical Recipes [22] convension where the factor is \(2/n\) for real FFTs. This is different from the datacondAPI where the normalization constant 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 [11] .
  4. All of these routines leave the input array undamaged. (Except for LALREAL4VectorFFT().)
  5. LALMalloc() is used by all the fftw routines.

Prototypes

REAL4FFTPlan * XLALCreateREAL4FFTPlan (UINT4 size, int fwdflg, int measurelvl)
 Returns a new REAL4FFTPlan A REAL4FFTPlan is required to perform a FFT that involves real data. More...
 
REAL4FFTPlan * XLALCreateForwardREAL4FFTPlan (UINT4 size, int measurelvl)
 Returns a new REAL4FFTPlan for a forward transform. More...
 
REAL4FFTPlan * XLALCreateReverseREAL4FFTPlan (UINT4 size, int measurelvl)
 Returns a new REAL4FFTPlan for a reverse transform. More...
 
void XLALDestroyREAL4FFTPlan (REAL4FFTPlan *plan)
 Destroys a REAL4FFTPlan. More...
 
int XLALREAL4ForwardFFT (COMPLEX8Vector *output, const REAL4Vector *input, const REAL4FFTPlan *plan)
 Performs a forward FFT of REAL4 data. More...
 
int XLALREAL4ReverseFFT (REAL4Vector *output, const COMPLEX8Vector *input, const REAL4FFTPlan *plan)
 Performs a reverse FFT of REAL4 data. More...
 
int XLALREAL4VectorFFT (REAL4Vector *_LAL_RESTRICT_ output, const REAL4Vector *_LAL_RESTRICT_ input, const REAL4FFTPlan *plan)
 Perform a REAL4Vector to REAL4Vector FFT. More...
 
int XLALREAL4PowerSpectrum (REAL4Vector *_LAL_RESTRICT_ spec, const REAL4Vector *_LAL_RESTRICT_ data, const REAL4FFTPlan *plan)
 Computes the power spectrum of REAL4 data. More...
 
REAL8FFTPlan * XLALCreateREAL8FFTPlan (UINT4 size, int fwdflg, int measurelvl)
 Returns a new REAL8FFTPlan. More...
 
REAL8FFTPlan * XLALCreateForwardREAL8FFTPlan (UINT4 size, int measurelvl)
 Returns a new REAL8FFTPlan for a forward transform. More...
 
REAL8FFTPlan * XLALCreateReverseREAL8FFTPlan (UINT4 size, int measurelvl)
 Returns a new REAL8FFTPlan for a reverse transform. More...
 
void XLALDestroyREAL8FFTPlan (REAL8FFTPlan *plan)
 Destroys a REAL8FFTPlan. More...
 
int XLALREAL8ForwardFFT (COMPLEX16Vector *output, const REAL8Vector *input, const REAL8FFTPlan *plan)
 Performs a forward FFT of REAL8 data. More...
 
int XLALREAL8ReverseFFT (REAL8Vector *output, const COMPLEX16Vector *input, const REAL8FFTPlan *plan)
 Performs a reverse FFT of REAL8 data. More...
 
int XLALREAL8VectorFFT (REAL8Vector *output, const REAL8Vector *input, const REAL8FFTPlan *plan)
 Perform a REAL8Vector to REAL8Vector FFT. More...
 
int XLALREAL8PowerSpectrum (REAL8Vector *spec, const REAL8Vector *data, const REAL8FFTPlan *plan)
 Computes the power spectrum of REAL8 data. More...
 

Data Structures

struct  REAL4FFTPlan
 Plan to perform FFT of REAL4 data. More...
 
struct  REAL8FFTPlan
 Plan to perform FFT of REAL8 data. More...
 

Macros

#define SINGLE_PRECISION
 
#define tagRealFFTPlan   tagREAL4FFTPlan
 
#define RealFFTPlan   REAL4FFTPlan
 

Files

file  RealFFTTest.c
 Tests the routines in RealFFT.h.
 

Error Codes

#define REALFFTH_ENULL   1
 Null pointer. More...
 
#define REALFFTH_ENNUL   2
 Non-null pointer. More...
 
#define REALFFTH_ESIZE   4
 Invalid input size. More...
 
#define REALFFTH_ESZMM   8
 Size mismatch. More...
 
#define REALFFTH_ESLEN   16
 Invalid/mismatched sequence lengths. More...
 
#define REALFFTH_ESAME   32
 Input/Output data vectors are the same. More...
 
#define REALFFTH_ESIGN   64
 Incorrect plan sign. More...
 
#define REALFFTH_EDATA   128
 Bad input data: DC/Nyquist should be real. More...
 
#define REALFFTH_EALOC   256
 Memory allocation failed. More...
 
#define REALFFTH_EFFTW   512
 Error in FFTW. More...
 
#define REALFFTH_ESNGL   1024
 FFTW library is not single-precision. More...
 
#define REALFFTH_EINTL   2048
 Error in Intel FFT library. More...
 

Function Documentation

◆ XLALCreateREAL4FFTPlan()

REAL4FFTPlan* XLALCreateREAL4FFTPlan ( UINT4  size,
int  fwdflg,
int  measurelvl 
)

Returns a new REAL4FFTPlan A REAL4FFTPlan is required to perform a FFT that involves real data.

A different plan is required for each size of the real data vector 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}\,x[j]\]

where N, the size of the transform, is the length of the vector x. A reverse transform performs

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

where N, the size of the transform, is the length of the vector y.

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 real 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 REAL4FFTPlan structure is returned upon successful completion. Otherwise, a NULL pointer is returned and xlalErrno is set to indicate the error.
Errors:
The XLALCreateREAL4Plan() 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.

◆ XLALCreateForwardREAL4FFTPlan()

REAL4FFTPlan* XLALCreateForwardREAL4FFTPlan ( UINT4  size,
int  measurelvl 
)

Returns a new REAL4FFTPlan for a forward transform.

A REAL4FFTPlan is required to perform a FFT that involves real data. A different plan is required for each size of the real data vector. A forward transform performs

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

where N, the size of the transform, is the length of the vector x.

Parameters
[in]sizeThe number of points in the real 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 REAL4FFTPlan structure is returned upon successful completion. Otherwise, a NULL pointer is returned and xlalErrno is set to indicate the error.
Errors:
The XLALCreateForwardREAL4Plan() 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 120 of file CudaRealFFT.c.

◆ XLALCreateReverseREAL4FFTPlan()

REAL4FFTPlan* XLALCreateReverseREAL4FFTPlan ( UINT4  size,
int  measurelvl 
)

Returns a new REAL4FFTPlan for a reverse transform.

A REAL4FFTPlan is required to perform a FFT that involves real data. A reverse transform performs

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

where N, the size of the transform, is the length of the vector y.

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 real 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 REAL4FFTPlan structure is returned upon successful completion. Otherwise, a NULL pointer is returned and xlalErrno is set to indicate the error.
Errors:
The XLALCreateReverseREAL4Plan() 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 130 of file CudaRealFFT.c.

◆ XLALDestroyREAL4FFTPlan()

void XLALDestroyREAL4FFTPlan ( REAL4FFTPlan *  plan)

Destroys a REAL4FFTPlan.

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

Definition at line 140 of file CudaRealFFT.c.

◆ XLALREAL4ForwardFFT()

int XLALREAL4ForwardFFT ( COMPLEX8Vector output,
const REAL4Vector input,
const REAL4FFTPlan *  plan 
)

Performs a forward FFT of REAL4 data.

This routine performs the transformation:

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

where N, the size of the transform, is the length of the vector x.

Note
Due to the reality of the input data x, the following identity holds for the complex FFT data: \(z[N-k]=z^\ast[k]\). Therefore, the length of the output vector is equal to \(\lfloor N/2\rfloor + 1\) since the remaining "negative" frequency components can be obtained from the "positive" frequency components.
Parameters
[out]outputThe complex data vector z of length [N/2] + 1 that results from the transform
[in]inputThe real data vector x of length to be transformed
[in]planThe FFT plan to use for the transform
Returns
0 upon successful completion or non-zero upon failure.
Errors:
The XLALREAL4ForwardFFT() function shall fail if:
  • [XLAL_EFAULT] A NULL pointer is provided as one of the arguments.
  • [XLAL_EINVAL] A argument is invalid or the plan is for a reverse transform.
  • [XLAL_EBADLEN] The input vector, output vector, and plan size are incompatible.
  • [XLAL_ENOMEM] Insufficient storage space is available.

Definition at line 161 of file CudaRealFFT.c.

◆ XLALREAL4ReverseFFT()

int XLALREAL4ReverseFFT ( REAL4Vector output,
const COMPLEX8Vector input,
const REAL4FFTPlan *  plan 
)

Performs a reverse FFT of REAL4 data.

This routine performs the transformation:

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

where N, the size of the transform, is the length of the vector y.

Note
Due to the reality of the output data y, the following identity holds for the complex data: \(z[N-k]=z^\ast[k]\). Therefore, the length of the input vector is equal to \(\lfloor N/2\rfloor + 1\) since the remaining "negative" frequency components can be obtained from the "positive" frequency components.
Parameters
[out]outputThe real data vector y of length N that results from the transform
[in]inputThe complex data vector z of length [N/2] + 1 to be transformed
[in]planThe FFT plan to use for the transform
Returns
0 upon successful completion or non-zero upon failure.
Errors:
The XLALREAL4ForwardFFT() function shall fail if:
  • [XLAL_EFAULT] A NULL pointer is provided as one of the arguments.
  • [XLAL_EINVAL] A argument is invalid or the plan is for a reverse transform.
  • [XLAL_EBADLEN] The input vector, output vector, and plan size are incompatible.
  • [XLAL_ENOMEM] Insufficient storage space is available.
  • [XLAL_EDOM] Domain error if the DC component of the input data, z[0], is not purely real or if the length of the output vector N is even and the Nyquist component of the input data, z[N/2], is not purely real.

Definition at line 198 of file CudaRealFFT.c.

◆ XLALREAL4VectorFFT()

int XLALREAL4VectorFFT ( REAL4Vector *_LAL_RESTRICT_  output,
const REAL4Vector *_LAL_RESTRICT_  input,
const REAL4FFTPlan *  plan 
)

Perform a REAL4Vector to REAL4Vector FFT.

This routine computes

\[y[k]=\left\{\begin{array}{ll}\Re z[k]&0\le k\le\lfloor N/2\rfloor\\\Im z[N-k]&\lfloor N/2\rfloor<k<N\end{array}\right.\]

where

\[z[k] = \sum_{j=0}^{N-1} e^{\mp2\pi ijk/N}\,x[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 vector x.

Parameters
[out]outputThe real output data vector y of length N
[in]inputThe input real data vector x 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 XLALREAL4VectorFFT() function shall fail if:
  • [XLAL_EFAULT] A NULL pointer is provided as one of the arguments.
  • [XLAL_EINVAL] A argument is invalid or the plan is for a reverse transform.
  • [XLAL_EBADLEN] The input vector, output vector, and plan size are incompatible.
  • [XLAL_ENOMEM] Insufficient storage space is available.

Definition at line 233 of file CudaRealFFT.c.

◆ XLALREAL4PowerSpectrum()

int XLALREAL4PowerSpectrum ( REAL4Vector *_LAL_RESTRICT_  spec,
const REAL4Vector *_LAL_RESTRICT_  data,
const REAL4FFTPlan *  plan 
)

Computes the power spectrum of REAL4 data.

This routine computes

\[P[k]=\left\{\begin{array}{ll}|z[0]|^2&k=0\\2|z[k]|^2&1\leq \lfloor (N+1)/2\rfloor\\|z[N/2]|^2&k=N/2,\;\mbox{$N$ even}\end{array}\right.\]

where

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

and N is the length of the input vector x.

Parameters
[out]specThe real power spectrum P of length [N/2] + 1 of the data x
[in]dataThe input real data vector x of length N
[in]planThe FFT plan to use for the transform
Returns
0 upon successful completion or non-zero upon failure.
Errors:
The XLALREAL4PowerSpectrum() 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.

◆ XLALCreateREAL8FFTPlan()

REAL8FFTPlan* XLALCreateREAL8FFTPlan ( UINT4  size,
int  fwdflg,
int  measurelvl 
)

Returns a new REAL8FFTPlan.

A REAL8FFTPlan is required to perform a FFT that involves real data. A different plan is required for each size of the real data vector 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}\,x[j]\]

where N, the size of the transform, is the length of the vector x. A reverse transform performs

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

where N, the size of the transform, is the length of the vector y.

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 real 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 REAL8FFTPlan structure is returned upon successful completion. Otherwise, a NULL pointer is returned and xlalErrno is set to indicate the error.
Errors:
The XLALCreateREAL8Plan() 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 366 of file CudaRealFFT.c.

◆ XLALCreateForwardREAL8FFTPlan()

REAL8FFTPlan* XLALCreateForwardREAL8FFTPlan ( UINT4  size,
int  measurelvl 
)

Returns a new REAL8FFTPlan for a forward transform.

A REAL8FFTPlan is required to perform a FFT that involves real data. A different plan is required for each size of the real data vector. A forward transform performs

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

where N, the size of the transform, is the length of the vector x.

Parameters
[in]sizeThe number of points in the real 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 REAL8FFTPlan structure is returned upon successful completion. Otherwise, a NULL pointer is returned and xlalErrno is set to indicate the error.
Errors:
The XLALCreateForwardREAL8Plan() 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 434 of file CudaRealFFT.c.

◆ XLALCreateReverseREAL8FFTPlan()

REAL8FFTPlan* XLALCreateReverseREAL8FFTPlan ( UINT4  size,
int  measurelvl 
)

Returns a new REAL8FFTPlan for a reverse transform.

A REAL8FFTPlan is required to perform a FFT that involves real data. A reverse transform performs

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

where N, the size of the transform, is the length of the vector y.

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 real 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 REAL8FFTPlan structure is returned upon successful completion. Otherwise, a NULL pointer is returned and xlalErrno is set to indicate the error.
Errors:
The XLALCreateReverseREAL8Plan() 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 444 of file CudaRealFFT.c.

◆ XLALDestroyREAL8FFTPlan()

void XLALDestroyREAL8FFTPlan ( REAL8FFTPlan *  plan)

Destroys a REAL8FFTPlan.

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

Definition at line 454 of file CudaRealFFT.c.

◆ XLALREAL8ForwardFFT()

int XLALREAL8ForwardFFT ( COMPLEX16Vector output,
const REAL8Vector input,
const REAL8FFTPlan *  plan 
)

Performs a forward FFT of REAL8 data.

This routine performs the transformation:

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

where N, the size of the transform, is the length of the vector x.

Note
Due to the reality of the input data x, the following identity holds for the complex FFT data: \(z[N-k]=z^\ast[k]\). Therefore, the length of the output vector is equal to \(\lfloor N/2\rfloor + 1\) since the remaining "negative" frequency components can be obtained from the "positive" frequency components.
Parameters
[out]outputThe complex data vector z of length [N/2] + 1 that results from the transform
[in]inputThe real data vector x of length to be transformed
[in]planThe FFT plan to use for the transform
Returns
0 upon successful completion or non-zero upon failure.
Errors:
The XLALREAL8ForwardFFT() function shall fail if:
  • [XLAL_EFAULT] A NULL pointer is provided as one of the arguments.
  • [XLAL_EINVAL] A argument is invalid or the plan is for a reverse transform.
  • [XLAL_EBADLEN] The input vector, output vector, and plan size are incompatible.
  • [XLAL_ENOMEM] Insufficient storage space is available.

Definition at line 471 of file CudaRealFFT.c.

◆ XLALREAL8ReverseFFT()

int XLALREAL8ReverseFFT ( REAL8Vector output,
const COMPLEX16Vector input,
const REAL8FFTPlan *  plan 
)

Performs a reverse FFT of REAL8 data.

This routine performs the transformation:

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

where N, the size of the transform, is the length of the vector y.

Note
Due to the reality of the output data y, the following identity holds for the complex data: \(z[N-k]=z^\ast[k]\). Therefore, the length of the input vector is equal to \(\lfloor N/2\rfloor + 1\) since the remaining "negative" frequency components can be obtained from the "positive" frequency components.
Parameters
[out]outputThe real data vector y of length N that results from the transform
[in]inputThe complex data vector z of length [N/2] + 1 to be transformed
[in]planThe FFT plan to use for the transform
Returns
0 upon successful completion or non-zero upon failure.
Errors:
The XLALREAL8ForwardFFT() function shall fail if:
  • [XLAL_EFAULT] A NULL pointer is provided as one of the arguments.
  • [XLAL_EINVAL] A argument is invalid or the plan is for a reverse transform.
  • [XLAL_EBADLEN] The input vector, output vector, and plan size are incompatible.
  • [XLAL_ENOMEM] Insufficient storage space is available.
  • [XLAL_EDOM] Domain error if the DC component of the input data, z[0], is not purely real or if the length of the output vector N is even and the Nyquist component of the input data, z[N/2], is not purely real.

Definition at line 518 of file CudaRealFFT.c.

◆ XLALREAL8VectorFFT()

int XLALREAL8VectorFFT ( REAL8Vector output,
const REAL8Vector input,
const REAL8FFTPlan *  plan 
)

Perform a REAL8Vector to REAL8Vector FFT.

This routine computes

\[y[k]=\left\{\begin{array}{ll}\Re z[k]&0\le k\le\lfloor N/2\rfloor\\\Im z[N-k]&\lfloor N/2\rfloor<k<N\end{array}\right.\]

where

\[z[k] = \sum_{j=0}^{N-1} e^{\mp2\pi ijk/N}\,x[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 vector x.

Parameters
[out]outputThe real output data vector y of length N
[in]inputThe input real data vector x 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 XLALREAL8VectorFFT() 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.

◆ XLALREAL8PowerSpectrum()

int XLALREAL8PowerSpectrum ( REAL8Vector spec,
const REAL8Vector data,
const REAL8FFTPlan *  plan 
)

Computes the power spectrum of REAL8 data.

This routine computes

\[P[k]=\left\{\begin{array}{ll}|z[0]|^2 & k=0\\2|z[k]|^2 & 1\leq \lfloor (N+1)/2\rfloor\\ |z[N/2]|^2 & k=N/2,\;\mbox{$N$ even}\end{array}\right.\]

where

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

and N is the length of the input vector x.

Parameters
[out]specThe real power spectrum P of length [N/2] + 1 of the data x
[in]dataThe input real data vector x of length N
[in]planThe FFT plan to use for the transform
Returns
0 upon successful completion or non-zero upon failure.
Errors:
The XLALREAL8PowerSpectrum() function shall fail if:
  • [XLAL_EFAULT] A NULL pointer is provided as one of the arguments.
  • [XLAL_EINVAL] A argument is invalid or the plan is for a reverse transform.
  • [XLAL_EBADLEN] The input vector, output vector, and plan size are incompatible.
  • [XLAL_ENOMEM] Insufficient storage space is available.

Definition at line 589 of file CudaRealFFT.c.

Macro Definition Documentation

◆ SINGLE_PRECISION

#define SINGLE_PRECISION

Definition at line 175 of file RealFFT.c.

◆ REALFFTH_ENULL

#define REALFFTH_ENULL   1

Null pointer.

Definition at line 161 of file RealFFT.h.

◆ REALFFTH_ENNUL

#define REALFFTH_ENNUL   2

Non-null pointer.

Definition at line 162 of file RealFFT.h.

◆ REALFFTH_ESIZE

#define REALFFTH_ESIZE   4

Invalid input size.

Definition at line 163 of file RealFFT.h.

◆ REALFFTH_ESZMM

#define REALFFTH_ESZMM   8

Size mismatch.

Definition at line 164 of file RealFFT.h.

◆ REALFFTH_ESLEN

#define REALFFTH_ESLEN   16

Invalid/mismatched sequence lengths.

Definition at line 165 of file RealFFT.h.

◆ REALFFTH_ESAME

#define REALFFTH_ESAME   32

Input/Output data vectors are the same.

Definition at line 166 of file RealFFT.h.

◆ REALFFTH_ESIGN

#define REALFFTH_ESIGN   64

Incorrect plan sign.

Definition at line 167 of file RealFFT.h.

◆ REALFFTH_EDATA

#define REALFFTH_EDATA   128

Bad input data: DC/Nyquist should be real.

Definition at line 168 of file RealFFT.h.

◆ REALFFTH_EALOC

#define REALFFTH_EALOC   256

Memory allocation failed.

Definition at line 169 of file RealFFT.h.

◆ REALFFTH_EFFTW

#define REALFFTH_EFFTW   512

Error in FFTW.

Definition at line 170 of file RealFFT.h.

◆ REALFFTH_ESNGL

#define REALFFTH_ESNGL   1024

FFTW library is not single-precision.

Definition at line 171 of file RealFFT.h.

◆ REALFFTH_EINTL

#define REALFFTH_EINTL   2048

Error in Intel FFT library.

Definition at line 172 of file RealFFT.h.

◆ tagRealFFTPlan

#define tagRealFFTPlan   tagREAL4FFTPlan

Definition at line 194 of file RealFFT.h.

◆ RealFFTPlan

#define RealFFTPlan   REAL4FFTPlan

Definition at line 195 of file RealFFT.h.