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

details Image Correlation VIGRA

Functions

template<... >
void crossCorrelation (...)
 This function performes a (slow) cross-correlation. More...
 
template<... >
void fastCrossCorrelation (...)
 This function performes a fast cross-correlation. More...
 
template<... >
void fastNormalizedCrossCorrelation (...)
 This function performes a fast normalized cross-correlation. More...
 
template<... >
void normalizedCrossCorrelation (...)
 This function performes a (slow) normalized cross-correlation. More...
 

Detailed Description

Fast (FFT-based) and slow algorithms to estimate the correlation between images.

Function Documentation

void vigra::fastCrossCorrelation (   ...)

This function performes a fast cross-correlation.

This function performes a fast cross-correlation using the Fast Fourier Transform and the dependency of the convolution and the correlation in Fourier space.

The input pixel type T1 must be a linear space over the window functions' value_type T, i.e. addition of source values, multiplication with functions' values, and NumericTraits must be defined. The mask's value_type must be an algebraic field, i.e. the arithmetic operations (+, -, *, /) and NumericTraits must be defined.

By default, the borders are filled with zeros. Use the clearBorders switch to change that behavior if you need to.

Declarations:

pass 2D array views:

namespace vigra {
template <class T1, class S1,
class T2, class S2,
class T3, class S3,
unsigned int N>
void
fastCrossCorrelation(MultiArrayView<N, T1, S1> const & in,
MultiArrayView<N, T2, S2> const & mask,
MultiArrayView<N, T3, S3> out,
bool clearBorders=true);
}

Usage:

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

unsigned int m_w=51, m_h=51;
unsigned int w=1000, h=1000;
MultiArray<2, float> mask(m_w,m_h), src(w,h), dest(w,h);
...
//compute fast cross correlation of mask and image -> dest
fastCrossCorrelation(mask, src, dest);

Preconditions:

The image must be larger than the size of the mask.

void vigra::fastNormalizedCrossCorrelation (   ...)

This function performes a fast normalized cross-correlation.

To compute the normalized cross correlation in a fast way, it is using the Fast Fourier Transform and sum-image look-up-tables as it is suggested by J.P.Lewis (1995): "Fast Normalized Cross-Correlation".

The input pixel type T1 must be a linear space over the window functions' value_type T, i.e. addition of source values, multiplication with functions' values, and NumericTraits must be defined. The mask's value_type must be an algebraic field, i.e. the arithmetic operations (+, -, *, /) and NumericTraits must be defined.

By default, the borders are filled with zeros. Use the clearBorders switch to change that behavior if you need to.

Declarations:

pass 2D array views:

namespace vigra {
template <class T1, class S1,
class T2, class S2,
class T3, class S3>
void
fastNormalizedCrossCorrelation(MultiArrayView<2, T1, S1> const & in,
MultiArrayView<2, T2, S2> const & mask,
MultiArrayView<2, T3, S3> out,
bool clearBorders=true);
}

Usage:

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

unsigned int m_w=51, m_h=51;
unsigned int w=1000, h=1000;
MultiArray<2, float> mask(m_w,m_h), src(w,h), dest(w,h);
...
//compute normalized cross correlation of mask and image -> dest
fastNormalizedCrossCorrelation(mask, src, dest);

Preconditions:

The image must be larger than the size of the mask.

void vigra::crossCorrelation (   ...)

This function performes a (slow) cross-correlation.

This function performes a (slow) cross-correlation using the window function environment and comparison of the mask with the underlying image part for each pixel. This may however be faster for very few comparisons.

The input pixel type T1 must be a linear space over the window functions' value_type T, i.e. addition of source values, multiplication with functions' values, and NumericTraits must be defined. The mask's value_type must be an algebraic field, i.e. the arithmetic operations (+, -, *, /) and NumericTraits must be defined.

Declarations:

pass 2D array views:

namespace vigra {
template <class T1, class S1,
class T2, class S2,
class T3, class S3>
void
crossCorrelation(MultiArrayView<2, T1, S1> const & in,
MultiArrayView<2, T2, S2> const & mask,
MultiArrayView<2, T3, S3> out);
}

show deprecated declarations

Usage:

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

unsigned int m_w=51, m_h=51;
unsigned int w=1000, h=1000;
MultiArray<2, float> mask(m_w,m_h), src(w,h), dest(w,h);
...
//compute (slow) cross correlation of mask and image -> dest
crossCorrelation(mask, src, dest);

Preconditions:

The image must be larger than the size of the mask.

void vigra::normalizedCrossCorrelation (   ...)

This function performes a (slow) normalized cross-correlation.

This function performes a (slow) normalized cross-correlation using the window function environment and comparison of the mask with the underlying image part for each pixel. This may however be faster for very few comparisons.

The input pixel type T1 must be a linear space over the window functions' value_type T, i.e. addition of source values, multiplication with functions' values, and NumericTraits must be defined. The mask's value_type must be an algebraic field, i.e. the arithmetic operations (+, -, *, /) and NumericTraits must be defined.

Declarations:

pass 2D array views:

namespace vigra {
template <class T1, class S1,
class T2, class S2,
class T3, class S3>
void
normalizedCrossCorrelation(MultiArrayView<2, T1, S1> const & in,
MultiArrayView<2, T2, S2> const & mask,
MultiArrayView<2, T3, S3> out);
}

show deprecated declarations

Usage:

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

unsigned int m_w=51, m_h=51;
unsigned int w=1000, h=1000;
MultiArray<2, float> mask(m_w,m_h), src(w,h), dest(w,h);
...
//compute (slow) normalized cross correlation of mask and image -> dest
normalizedCrossCorrelation(mask, src, dest);

Preconditions:

The image must be larger than the size of the mask.

© 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)