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

details FFTWConvolvePlan< N, Real > Class Template Reference VIGRA

#include <vigra/multi_fft.hxx>

Inheritance diagram for FFTWConvolvePlan< N, Real >:
FFTWCorrelatePlan< N, Real >

Public Member Functions

template<class C1 , class C2 , class C3 >
void execute (MultiArrayView< N, Real, C1 > in, MultiArrayView< N, Real, C2 > kernel, MultiArrayView< N, Real, C3 > out)
 Execute a plan to convolve a real array with a real kernel. More...
 
template<class C1 , class C2 , class C3 >
void execute (MultiArrayView< N, Real, C1 > in, MultiArrayView< N, FFTWComplex< Real >, C2 > kernel, MultiArrayView< N, Real, C3 > out)
 Execute a plan to convolve a real array with a complex kernel. More...
 
template<class C1 , class C2 , class C3 >
void execute (MultiArrayView< N, FFTWComplex< Real >, C1 > in, MultiArrayView< N, FFTWComplex< Real >, C2 > kernel, MultiArrayView< N, FFTWComplex< Real >, C3 > out)
 Execute a plan to convolve a complex array with a complex kernel. More...
 
template<class C1 , class KernelIterator , class OutIterator >
void executeMany (MultiArrayView< N, FFTWComplex< Real >, C1 > in, KernelIterator kernels, KernelIterator kernelsEnd, OutIterator outs)
 Execute a plan to convolve a complex array with a sequence of kernels. More...
 
template<class C1 , class KernelIterator , class OutIterator >
void executeMany (MultiArrayView< N, Real, C1 > in, KernelIterator kernels, KernelIterator kernelsEnd, OutIterator outs)
 Execute a plan to convolve a real array with a sequence of kernels. More...
 
 FFTWConvolvePlan ()
 Create an empty plan. More...
 
template<class C1 , class C2 , class C3 >
 FFTWConvolvePlan (MultiArrayView< N, Real, C1 > in, MultiArrayView< N, Real, C2 > kernel, MultiArrayView< N, Real, C3 > out, unsigned int planner_flags=FFTW_ESTIMATE)
 Create a plan to convolve a real array with a real kernel. More...
 
template<class C1 , class C2 , class C3 >
 FFTWConvolvePlan (MultiArrayView< N, Real, C1 > in, MultiArrayView< N, FFTWComplex< Real >, C2 > kernel, MultiArrayView< N, Real, C3 > out, unsigned int planner_flags=FFTW_ESTIMATE)
 Create a plan to convolve a real array with a complex kernel. More...
 
template<class C1 , class C2 , class C3 >
 FFTWConvolvePlan (MultiArrayView< N, FFTWComplex< Real >, C1 > in, MultiArrayView< N, FFTWComplex< Real >, C2 > kernel, MultiArrayView< N, FFTWComplex< Real >, C3 > out, bool fourierDomainKernel, unsigned int planner_flags=FFTW_ESTIMATE)
 Create a plan to convolve a complex array with a complex kernel. More...
 
template<class C1 , class C2 , class C3 >
 FFTWConvolvePlan (Shape inOut, Shape kernel, bool useFourierKernel=false, unsigned int planner_flags=FFTW_ESTIMATE)
 Create a plan from just the shape information. More...
 
template<class C1 , class C2 , class C3 >
void init (MultiArrayView< N, Real, C1 > in, MultiArrayView< N, Real, C2 > kernel, MultiArrayView< N, Real, C3 > out, unsigned int planner_flags=FFTW_ESTIMATE)
 Init a plan to convolve a real array with a real kernel. More...
 
template<class C1 , class C2 , class C3 >
void init (MultiArrayView< N, Real, C1 > in, MultiArrayView< N, FFTWComplex< Real >, C2 > kernel, MultiArrayView< N, Real, C3 > out, unsigned int planner_flags=FFTW_ESTIMATE)
 Init a plan to convolve a real array with a complex kernel. More...
 
template<class C1 , class C2 , class C3 >
void init (MultiArrayView< N, FFTWComplex< Real >, C1 > in, MultiArrayView< N, FFTWComplex< Real >, C2 > kernel, MultiArrayView< N, FFTWComplex< Real >, C3 > out, bool fourierDomainKernel, unsigned int planner_flags=FFTW_ESTIMATE)
 Init a plan to convolve a complex array with a complex kernel. More...
 
template<class C1 , class KernelIterator , class OutIterator >
void initMany (MultiArrayView< N, Real, C1 > in, KernelIterator kernels, KernelIterator kernelsEnd, OutIterator outs, unsigned int planner_flags=FFTW_ESTIMATE)
 Init a plan to convolve a real array with a sequence of kernels. More...
 
template<class C1 , class KernelIterator , class OutIterator >
void initMany (MultiArrayView< N, FFTWComplex< Real >, C1 > in, KernelIterator kernels, KernelIterator kernelsEnd, OutIterator outs, bool fourierDomainKernels, unsigned int planner_flags=FFTW_ESTIMATE)
 Init a plan to convolve a complex array with a sequence of kernels. More...
 

Detailed Description

template<unsigned int N, class Real>
class vigra::FFTWConvolvePlan< N, Real >

C++ wrapper for a pair of FFTW plans used to perform FFT-based convolution.

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. It always creates a pair of plans, one for the forward and one for the inverse transform required for convolution.

Usually, you use this class only indirectly via convolveFFT() and its variants. 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 convolutions.

Usage:

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

// convolve a real array with a real kernel
MultiArray<2, double> src(Shape2(w, h)), dest(Shape2(w, h));
MultiArray<2, double> spatial_kernel(Shape2(9, 9));
Gaussian<double> gauss(1.0);
for(int y=0; y<9; ++y)
for(int x=0; x<9; ++x)
spatial_kernel(x, y) = gauss(x-4.0)*gauss(y-4.0);
// create an optimized plan by measuring the speed of several algorithm variants
FFTWConvolvePlan<2, double> plan(src, spatial_kernel, dest, FFTW_MEASURE);
plan.execute(src, spatial_kernel, dest);

Constructor & Destructor Documentation

Create an empty plan.

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

FFTWConvolvePlan ( MultiArrayView< N, Real, C1 >  in,
MultiArrayView< N, Real, C2 >  kernel,
MultiArrayView< N, Real, C3 >  out,
unsigned int  planner_flags = FFTW_ESTIMATE 
)

Create a plan to convolve a real array with a real kernel.

The kernel must be defined in the spatial domain. See convolveFFT() for detailed information on required shapes and internal padding.

  • 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".
FFTWConvolvePlan ( MultiArrayView< N, Real, C1 >  in,
MultiArrayView< N, FFTWComplex< Real >, C2 >  kernel,
MultiArrayView< N, Real, C3 >  out,
unsigned int  planner_flags = FFTW_ESTIMATE 
)

Create a plan to convolve a real array with a complex kernel.

The kernel must be defined in the Fourier domain, using the half-space format. See convolveFFT() for detailed information on required shapes and internal padding.

  • 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".
FFTWConvolvePlan ( MultiArrayView< N, FFTWComplex< Real >, C1 >  in,
MultiArrayView< N, FFTWComplex< Real >, C2 >  kernel,
MultiArrayView< N, FFTWComplex< Real >, C3 >  out,
bool  fourierDomainKernel,
unsigned int  planner_flags = FFTW_ESTIMATE 
)

Create a plan to convolve a complex array with a complex kernel.

See convolveFFT() for detailed information on required shapes and internal padding.

  • fourierDomainKernel determines if the kernel is defined in the spatial or Fourier domain.
  • 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".
FFTWConvolvePlan ( Shape  inOut,
Shape  kernel,
bool  useFourierKernel = false,
unsigned int  planner_flags = FFTW_ESTIMATE 
)

Create a plan from just the shape information.

See convolveFFT() for detailed information on required shapes and internal padding.

  • fourierDomainKernel determines if the kernel is defined in the spatial or Fourier domain.
  • 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, Real, C1 >  in,
MultiArrayView< N, Real, C2 >  kernel,
MultiArrayView< N, Real, C3 >  out,
unsigned int  planner_flags = FFTW_ESTIMATE 
)

Init a plan to convolve a real array with a real kernel.

See the constructor with the same signature for details.

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

Init a plan to convolve a real array with a complex kernel.

See the constructor with the same signature for details.

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

Init a plan to convolve a complex array with a complex kernel.

See the constructor with the same signature for details.

void initMany ( MultiArrayView< N, Real, C1 >  in,
KernelIterator  kernels,
KernelIterator  kernelsEnd,
OutIterator  outs,
unsigned int  planner_flags = FFTW_ESTIMATE 
)

Init a plan to convolve a real array with a sequence of kernels.

The kernels can be either real or complex. The sequences kernels and outs must have the same length. See the corresponding constructors for single kernels for details.

void initMany ( MultiArrayView< N, FFTWComplex< Real >, C1 >  in,
KernelIterator  kernels,
KernelIterator  kernelsEnd,
OutIterator  outs,
bool  fourierDomainKernels,
unsigned int  planner_flags = FFTW_ESTIMATE 
)

Init a plan to convolve a complex array with a sequence of kernels.

The kernels must be complex as well. The sequences kernels and outs must have the same length. See the corresponding constructors for single kernels for details.

void execute ( MultiArrayView< N, Real, C1 >  in,
MultiArrayView< N, Real, C2 >  kernel,
MultiArrayView< N, Real, C3 >  out 
)

Execute a plan to convolve a real array with a real kernel.

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 >  kernel,
MultiArrayView< N, Real, C3 >  out 
)

Execute a plan to convolve a real array with a complex kernel.

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, FFTWComplex< Real >, C2 >  kernel,
MultiArrayView< N, FFTWComplex< Real >, C3 >  out 
)

Execute a plan to convolve a complex array with a complex kernel.

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 executeMany ( MultiArrayView< N, FFTWComplex< Real >, C1 >  in,
KernelIterator  kernels,
KernelIterator  kernelsEnd,
OutIterator  outs 
)

Execute a plan to convolve a complex array with a sequence of kernels.

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

void executeMany ( MultiArrayView< N, Real, C1 >  in,
KernelIterator  kernels,
KernelIterator  kernelsEnd,
OutIterator  outs 
)

Execute a plan to convolve a real array with a sequence of kernels.

The array shapes must be the same as in the corresponding init function or constructor. However, executeMany() 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)