[ VIGRA Homepage | Function Index | Class Index | Namespaces | File List | Main Page ]

details FFTWPlan< N, Real > Class Template Reference VIGRA

#include <vigra/multi_fft.hxx>

Public Member Functions

template<class C1 , class C2 >
void execute (MultiArrayView< N, FFTWComplex< Real >, C1 > in, MultiArrayView< N, FFTWComplex< Real >, C2 > out) const
 Execute a complex-to-complex transform. More...
 
template<class C1 , class C2 >
void execute (MultiArrayView< N, Real, C1 > in, MultiArrayView< N, FFTWComplex< Real >, C2 > out) const
 Execute a real-to-complex transform. More...
 
template<class C1 , class C2 >
void execute (MultiArrayView< N, FFTWComplex< Real >, C1 > in, MultiArrayView< N, Real, C2 > out) const
 Execute a complex-to-real transform. More...
 
 FFTWPlan ()
 Create an empty plan. More...
 
template<class C1 , class C2 >
 FFTWPlan (MultiArrayView< N, FFTWComplex< Real >, C1 > in, MultiArrayView< N, FFTWComplex< Real >, C2 > out, int SIGN, unsigned int planner_flags=FFTW_ESTIMATE)
 Create a plan for a complex-to-complex transform. More...
 
template<class C1 , class C2 >
 FFTWPlan (MultiArrayView< N, Real, C1 > in, MultiArrayView< N, FFTWComplex< Real >, C2 > out, unsigned int planner_flags=FFTW_ESTIMATE)
 Create a plan for a real-to-complex transform. More...
 
template<class C1 , class C2 >
 FFTWPlan (MultiArrayView< N, FFTWComplex< Real >, C1 > in, MultiArrayView< N, Real, C2 > out, unsigned int planner_flags=FFTW_ESTIMATE)
 Create a plan for a complex-to-real transform. More...
 
 FFTWPlan (FFTWPlan const &other)
 Copy constructor.
 
template<class C1 , class C2 >
void init (MultiArrayView< N, FFTWComplex< Real >, C1 > in, MultiArrayView< N, FFTWComplex< Real >, C2 > out, int SIGN, unsigned int planner_flags=FFTW_ESTIMATE)
 Init a complex-to-complex transform. More...
 
template<class C1 , class C2 >
void init (MultiArrayView< N, Real, C1 > in, MultiArrayView< N, FFTWComplex< Real >, C2 > out, unsigned int planner_flags=FFTW_ESTIMATE)
 Init a real-to-complex transform. More...
 
template<class C1 , class C2 >
void init (MultiArrayView< N, FFTWComplex< Real >, C1 > in, MultiArrayView< N, Real, C2 > out, unsigned int planner_flags=FFTW_ESTIMATE)
 Init a complex-to-real transform. More...
 
FFTWPlanoperator= (FFTWPlan const &other)
 Copy assigment.
 
 ~FFTWPlan ()
 Destructor.
 

Detailed Description

template<unsigned int N, class Real = double>
class vigra::FFTWPlan< N, Real >

C++ wrapper for FFTW plans.

The class encapsulates the calls to fftw_plan_dft_2d, fftw_execute, and fftw_destroy_plan (and their float and long double counterparts) in an easy-to-use interface.

Usually, you use this class only indirectly via fourierTransform() and fourierTransformInverse(). You only need this class if you want to have more control about FFTW's planning process (by providing non-default planning flags) and/or want to re-use plans for several transformations.

Usage:

#include <vigra/multi_fft.hxx>
Namespace: vigra

// compute complex Fourier transform of a real image
MultiArray<2, double> src(Shape2(w, h));
MultiArray<2, FFTWComplex<double> > fourier(Shape2(w, h));
// create an optimized plan by measuring the speed of several algorithm variants
FFTWPlan<2, double> plan(src, fourier, FFTW_MEASURE);
plan.execute(src, fourier);

Constructor & Destructor Documentation

FFTWPlan ( )

Create an empty plan.

The plan can be initialized later by one of the init() functions.

FFTWPlan ( MultiArrayView< N, FFTWComplex< Real >, C1 >  in,
MultiArrayView< N, FFTWComplex< Real >, C2 >  out,
int  SIGN,
unsigned int  planner_flags = FFTW_ESTIMATE 
)

Create a plan for a complex-to-complex transform.

  • SIGN must be FFTW_FORWARD or FFTW_BACKWARD according to the desired transformation direction.
  • planner_flags must be a combination of the planner flags defined by the FFTW library. The default FFTW_ESTIMATE will guess optimal algorithm settings or read them from pre-loaded "wisdom".
FFTWPlan ( MultiArrayView< N, Real, C1 >  in,
MultiArrayView< N, FFTWComplex< Real >, C2 >  out,
unsigned int  planner_flags = FFTW_ESTIMATE 
)

Create a plan for a real-to-complex transform.

This always refers to a forward transform. The shape of the output determines if a standard transform (when out.shape() == in.shape()) or an R2C transform (when out.shape() == fftwCorrespondingShapeR2C(in.shape())) will be executed.

  • planner_flags must be a combination of the planner flags defined by the FFTW library. The default FFTW_ESTIMATE will guess optimal algorithm settings or read them from pre-loaded "wisdom".
FFTWPlan ( MultiArrayView< N, FFTWComplex< Real >, C1 >  in,
MultiArrayView< N, Real, C2 >  out,
unsigned int  planner_flags = FFTW_ESTIMATE 
)

Create a plan for a complex-to-real transform.

This always refers to a inverse transform. The shape of the input determines if a standard transform (when in.shape() == out.shape()) or a C2R transform (when in.shape() == fftwCorrespondingShapeR2C(out.shape())) will be executed.

  • planner_flags must be a combination of the planner flags defined by the FFTW library. The default FFTW_ESTIMATE will guess optimal algorithm settings or read them from pre-loaded "wisdom".

Member Function Documentation

void init ( MultiArrayView< N, FFTWComplex< Real >, C1 >  in,
MultiArrayView< N, FFTWComplex< Real >, C2 >  out,
int  SIGN,
unsigned int  planner_flags = FFTW_ESTIMATE 
)

Init a complex-to-complex transform.

See the constructor with the same signature for details.

void init ( MultiArrayView< N, Real, C1 >  in,
MultiArrayView< N, FFTWComplex< Real >, C2 >  out,
unsigned int  planner_flags = FFTW_ESTIMATE 
)

Init a real-to-complex transform.

See the constructor with the same signature for details.

void init ( MultiArrayView< N, FFTWComplex< Real >, C1 >  in,
MultiArrayView< N, Real, C2 >  out,
unsigned int  planner_flags = FFTW_ESTIMATE 
)

Init a complex-to-real transform.

See the constructor with the same signature for details.

void execute ( MultiArrayView< N, FFTWComplex< Real >, C1 >  in,
MultiArrayView< N, FFTWComplex< Real >, C2 >  out 
) const

Execute a complex-to-complex transform.

The array shapes must be the same as in the corresponding init function or constructor. However, execute() can be called several times on the same plan, even with different arrays, as long as they have the appropriate shapes.

void execute ( MultiArrayView< N, Real, C1 >  in,
MultiArrayView< N, FFTWComplex< Real >, C2 >  out 
) const

Execute a real-to-complex transform.

The array shapes must be the same as in the corresponding init function or constructor. However, execute() can be called several times on the same plan, even with different arrays, as long as they have the appropriate shapes.

void execute ( MultiArrayView< N, FFTWComplex< Real >, C1 >  in,
MultiArrayView< N, Real, C2 >  out 
) const

Execute a complex-to-real transform.

The array shapes must be the same as in the corresponding init function or constructor. However, execute() can be called several times on the same plan, even with different arrays, as long as they have the appropriate shapes.


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

© Ullrich Köthe (ullrich.koethe@iwr.uni-heidelberg.de)
Heidelberg Collaboratory for Image Processing, University of Heidelberg, Germany

html generated using doxygen and Python
vigra 1.11.1 (Fri May 19 2017)