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

details FFTWCorrelatePlan< N, Real > Class Template Reference VIGRA

#include <vigra/multi_fft.hxx>

Inheritance diagram for FFTWCorrelatePlan< N, Real >:
FFTWConvolvePlan< 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 correlate a real array with a real kernel. More...
 
 FFTWCorrelatePlan ()
 Create an empty plan. More...
 
template<class C1 , class C2 , class C3 >
 FFTWCorrelatePlan (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 correlate a real array with a real kernel. More...
 
template<class C1 , class C2 , class C3 >
 FFTWCorrelatePlan (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...
 

Detailed Description

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

Like FFTWConvolvePlan, but performs correlation rather than convolution.

See vigra::FFTWConvolvePlan for details.

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
FFTWCorrelatePlan<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.

FFTWCorrelatePlan ( 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 correlate a real array with a real kernel.

The kernel must be defined in the spatial domain. See correlateFFT() 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".
FFTWCorrelatePlan ( 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 execute ( MultiArrayView< N, Real, C1 >  in,
MultiArrayView< N, Real, C2 >  kernel,
MultiArrayView< N, Real, C3 >  out 
)

Execute a plan to correlate 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.


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)