This section describes the API of castor-fftw.

Main API

fft

matrix<std::complex<float>> castor::fftw::fft(matrix<float> &X, int dim = 1)
matrix<std::complex<float>> castor::fftw::fft(matrix<std::complex<float>> &X, int dim = 1)
matrix<std::complex<double>> castor::fftw::fft(matrix<double> &X, int dim = 1)
matrix<std::complex<double>> castor::fftw::fft(matrix<std::complex<double>> &X, int dim = 1)

Computes the forward non-normalized 1d-Fourier transform. The default fft(A) assumes that the input A is one-dimensional (either a line or a column), even if it is not the case (in that case, the flattened matrix is taken into account).

std::size_t N=10;
matrix<double> A=rand<>(1,N);
matrix<std::complex<double>> Ahat = fft(A)/static_cast<double>(N);
It is also possible to compute the transform along the columns or the lines
std::size_t N=10, M=5;
matrix<double> A=rand<>(M,N);
auto Ahat_col = fft(A)/static_cast<double>(M); // ... = fft(A,1) !
auto Ahat_lin = fft(A,2)/static_cast<double>(N);

See also ifft.

ifft

matrix<std::complex<float>> castor::fftw::ifft(matrix<float> &X, int dim = 1)
matrix<std::complex<float>> castor::fftw::ifft(matrix<std::complex<float>> &X, int dim = 1)
matrix<std::complex<double>> castor::fftw::ifft(matrix<double> &X, int dim = 1)
matrix<std::complex<double>> castor::fftw::ifft(matrix<std::complex<double>> &X, int dim = 1)

Computes the backward 1d-Fourier transform. The interface is similar to fft().

See also fft.

fft2

matrix<std::complex<float>> castor::fftw::fft2(std::size_t m, std::size_t n, matrix<float> &X)
matrix<std::complex<float>> castor::fftw::fft2(std::size_t m, std::size_t n, matrix<std::complex<float>> &X)
matrix<std::complex<float>> castor::fftw::fft2(matrix<float> &X)
matrix<std::complex<float>> castor::fftw::fft2(matrix<std::complex<float>> &X)
matrix<std::complex<double>> castor::fftw::fft2(std::size_t m, std::size_t n, matrix<double> &X)
matrix<std::complex<double>> castor::fftw::fft2(std::size_t m, std::size_t n, matrix<std::complex<double>> &X)
matrix<std::complex<double>> castor::fftw::fft2(matrix<double> &X)
matrix<std::complex<double>> castor::fftw::fft2(matrix<std::complex<double>> &X)

Computes the forward non-normalized 2d-Fourier tranform.

std::size_t M=10, N=20;
matrix<std::complex<double>> A = rand<>(M,N) + M_1I*rand<>(M,N);
auto Ahat = fft2(A)/static_cast<double>(M*N);
It is also possible to give the input as a 1d or 2d array and provide its ‘true’ dimensions. The example above is equivalent to
auto Ahat = fft2(M,N,A)/static_cast<double>(M*N);

See also ifft2.

ifft2

matrix<std::complex<float>> castor::fftw::ifft2(std::size_t m, std::size_t n, matrix<float> &X)
matrix<std::complex<float>> castor::fftw::ifft2(std::size_t m, std::size_t n, matrix<std::complex<float>> &X)
matrix<std::complex<float>> castor::fftw::ifft2(matrix<float> &X)
matrix<std::complex<float>> castor::fftw::ifft2(matrix<std::complex<float>> &X)
matrix<std::complex<double>> castor::fftw::ifft2(std::size_t m, std::size_t n, matrix<double> &X)
matrix<std::complex<double>> castor::fftw::ifft2(std::size_t m, std::size_t n, matrix<std::complex<double>> &X)
matrix<std::complex<double>> castor::fftw::ifft2(matrix<double> &X)
matrix<std::complex<double>> castor::fftw::ifft2(matrix<std::complex<double>> &X)

Compute the backward 2d-Fourier transform.

See also fft2.

Helper functions

fftfreq

template<typename T>
matrix<T> castor::fftw::fftfreq(std::size_t n, T d)

Returns the Discrete Fourier Transform sample frequencies.

float       dt = 0.1;
std::size_t n  = 100;
//
matrix<float> freq = fftfreq<float>(n,dt);

fftshift

template<typename T>
matrix<T> castor::fftw::fftshift(matrix<T> const &A, int dim = 0)

Shifts the zero frequency component to the center of the spectrum.

matrix<float> freq = fftfreq(10,0.1);
auto freq_shifted = fftshift(freq);
The dim argument enables shifting along a given dimension.

See ifftshift, fft, fft2.

ifftshift

template<typename T>
matrix<T> castor::fftw::ifftshift(matrix<T> const &As, int dim = 0)

Inverse transform of fftshift().

See fftshift, ifft, ifft2.