Lattice Builder Manual
Software Package for Constructing Rank-1 Lattices
fftw< T > Struct Template Reference

Wrapper for a subset of FFTW: FFT's for real functions in one dimension. More...

#include <::c_api>

Classes

class  allocator
 STL allocator replacement using FFTW's memory allocation functions. More...
 
struct  c_api
 Specialization of c_api for double precision. More...
 

Public Types

typedef T real
 Real number.
 
typedef std::complex< T > complex
 Complex number.
 
typedef std::vector< real, allocator< real > > real_vector
 Vector of real numbers using FFTW allocation.
 
typedef std::vector< complex, allocator< complex > > complex_vector
 Vector of complex numbers using FFTW allocation.
 

Static Public Member Functions

static complex_vectorfft (const real_vector &v, complex_vector &result)
 Computes the real-to-complex Fourier transform of v into result. More...
 
static complex_vector::size_type fft_size (const real_vector &v)
 Returns the output size (exact) for real-to-complex Fourier transform of v. More...
 
static complex_vector fft (const real_vector &v)
 Computes the real-to-complex Fourier transform of v. More...
 
static real_vectorifft (const complex_vector &v, real_vector &result, bool normalize=true)
 Computes the complex-to-real Fourier transform of v into result. More...
 
static real_vector::size_type real_ifft_size (const complex_vector &v)
 Returns the output size (guessed) for complex-to-real Fourier transform of v. More...
 
static real_vector real_ifft (const complex_vector &v, bool normalize=true)
 Computes the complex-to-real Fourier transform of v. More...
 

Detailed Description

template<typename T>
struct fftw< T >

Wrapper for a subset of FFTW: FFT's for real functions in one dimension.

Specialization of for float precision.

Member Function Documentation

template<typename T>
static complex_vector& fftw< T >::fft ( const real_vector v,
complex_vector result 
)
inlinestatic

Computes the real-to-complex Fourier transform of v into result.

The size of the transform is that of the real component. The expected size of result is v.size() / 2 + 1. Because result is Hermitian, only half of the elements (plus one) are stored.

References fftw< T >::fft_size().

Referenced by fftw< T >::fft().

template<typename T>
static complex_vector fftw< T >::fft ( const real_vector v)
inlinestatic

Computes the real-to-complex Fourier transform of v.

See also
transform(const real_vector&, complex_vector&)

References fftw< T >::fft(), and fftw< T >::fft_size().

template<typename T>
static complex_vector::size_type fftw< T >::fft_size ( const real_vector v)
inlinestatic

Returns the output size (exact) for real-to-complex Fourier transform of v.

See also
transform(const real_vector&, complex_vector&)

Referenced by fftw< T >::fft(), and fftw< T >::ifft().

template<typename T>
static real_vector& fftw< T >::ifft ( const complex_vector v,
real_vector result,
bool  normalize = true 
)
inlinestatic

Computes the complex-to-real Fourier transform of v into result.

The size of the transform is that of the real component. The expected size of v is is result.size() / 2 + 1. Because v is Hermitian, only half of the elements (plus one) are considered.

Warning
If normalize is false, the result vector is not normalized; its components must be divided by result.size() for proper normalization.

References fftw< T >::fft_size().

Referenced by fftw< T >::real_ifft(), and LatBuilder::Kernel::RAlpha::valuesVector().

template<typename T>
static real_vector fftw< T >::real_ifft ( const complex_vector v,
bool  normalize = true 
)
inlinestatic

Computes the complex-to-real Fourier transform of v.

Assumes the size of the transform is 2 * (v.size() - 1). Returns a vector of size 2 * (v.size() - 1).

See also
transform(const complex_vector&, real_vector&)

References fftw< T >::ifft(), and fftw< T >::real_ifft_size().

template<typename T>
static real_vector::size_type fftw< T >::real_ifft_size ( const complex_vector v)
inlinestatic

Returns the output size (guessed) for complex-to-real Fourier transform of v.

See also
transform(const complex_vector&, real_vector&)

Referenced by fftw< T >::real_ifft().


The documentation for this struct was generated from the following file: