[ VIGRA Homepage | Function Index | Class Index | Namespaces | File List | Main Page ]
FFTWPlan< N, Real > Class Template Reference |
#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... | |
FFTWPlan & | operator= (FFTWPlan const &other) |
Copy assigment. | |
~FFTWPlan () | |
Destructor. | |
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
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.
FFTW_FORWARD
or FFTW_BACKWARD
according to the desired transformation direction. 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.
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.
FFTW_ESTIMATE
will guess optimal algorithm settings or read them from pre-loaded "wisdom". 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.
© Ullrich Köthe (ullrich.koethe@iwr.uni-heidelberg.de) |
html generated using doxygen and Python
|