[ VIGRA Homepage | Class Index | Function Index | File Index | Main Page ]
|
Common Filters | ![]() |
|---|
Functions | |
| template<...> void | convolveImage (SrcIterator supperleft, SrcIterator slowerright, SrcAccessor sa, DestIterator dupperleft, DestAccessor da, Kernel1D< T > const &kx, Kernel1D< T > const &ky) |
| Apply two separable filters successively, the first in x-direction, the second in y-direction. | |
| template<...> void | simpleSharpening (SrcIterator src_ul, SrcIterator src_lr, SrcAccessor src_acc, DestIterator dest_ul, DestAccessor dest_acc, double sharpening_factor) |
| Perform simple sharpening function. | |
| template<...> void | gaussianSharpening (SrcIterator src_ul, SrcIterator src_lr, SrcAccessor src_acc, DestIterator dest_ul, DestAccessor dest_acc, double sharpening_factor, double scale) |
| Perform sharpening function with gaussian filter. | |
| template<...> void | gaussianSmoothing (SrcIterator supperleft, SrcIterator slowerright, SrcAccessor sa, DestIterator dupperleft, DestAccessor da, double scale) |
| Perform isotropic Gaussian convolution. | |
| template<...> void | gaussianGradient (SrcIterator supperleft, SrcIterator slowerright, SrcAccessor sa, DestIteratorX dupperleftx, DestAccessorX dax, DestIteratorY dupperlefty, DestAccessorY day, double scale) |
| Calculate the gradient vector by means of a 1st derivatives of Gaussian filter. | |
| template<...> void | laplacianOfGaussian (SrcIterator supperleft, SrcIterator slowerright, SrcAccessor sa, DestIterator dupperleft, DestAccessor da, double scale) |
| Filter image with the Laplacian of Gaussian operator at the given scale. | |
| template<...> void | hessianMatrixOfGaussian (SrcIterator supperleft, SrcIterator slowerright, SrcAccessor sa, DestIteratorX dupperleftx, DestAccessorX dax, DestIteratorXY dupperleftxy, DestAccessorXY daxy, DestIteratorY dupperlefty, DestAccessorY day, double scale) |
| Filter image with the 2nd derivatives of the Gaussian at the given scale to get the Hessian matrix. | |
| template<...> void | structureTensor (SrcIterator supperleft, SrcIterator slowerright, SrcAccessor sa, DestIteratorX dupperleftx, DestAccessorX dax, DestIteratorXY dupperleftxy, DestAccessorXY daxy, DestIteratorY dupperlefty, DestAccessorY day, double inner_scale, double outer_scale) |
| Calculate the Structure Tensor for each pixel of and image, using Gaussian (derivative) filters. | |
Detailed Description |
|
void convolveImage (...) |
|
Apply two separable filters successively, the first in x-direction, the second in y-direction. This function is a shorthand for the concatenation of a call to separableConvolveX() and separableConvolveY() with the given kernels. Declarations: pass arguments explicitly: namespace vigra { template <class SrcIterator, class SrcAccessor, class DestIterator, class DestAccessor, class T> void convolveImage(SrcIterator supperleft, SrcIterator slowerright, SrcAccessor sa, DestIterator dupperleft, DestAccessor da, Kernel1D<T> const & kx, Kernel1D<T> const & ky); } use argument objects in conjuction with Argument Object Factories: namespace vigra { template <class SrcIterator, class SrcAccessor, class DestIterator, class DestAccessor, class T> inline void convolveImage(triple<SrcIterator, SrcIterator, SrcAccessor> src, pair<DestIterator, DestAccessor> dest, Kernel1D<T> const & kx, Kernel1D<T> const & ky); } Usage: #include "vigra/convolution.hxx"
vigra::FImage src(w,h), dest(w,h);
...
// implement sobel filter in x-direction
Kernel1D<double> kx, ky;
kx.initSymmetricGradient();
ky.initBinomial(1);
vigra::convolveImage(srcImageRange(src), destImage(dest), kx, ky);
|
|
void gaussianGradient (...) |
|
Calculate the gradient vector by means of a 1st derivatives of Gaussian filter. This function is a shorthand for the concatenation of a call to separableConvolveX() and separableConvolveY() with the appropriate kernels at the given scale. Note that this function produces two result images. Declarations: pass arguments explicitly: namespace vigra { template <class SrcIterator, class SrcAccessor, class DestIteratorX, class DestAccessorX, class DestIteratorY, class DestAccessorY> void gaussianGradient(SrcIterator supperleft, SrcIterator slowerright, SrcAccessor sa, DestIteratorX dupperleftx, DestAccessorX dax, DestIteratorY dupperlefty, DestAccessorY day, double scale); } use argument objects in conjuction with Argument Object Factories: namespace vigra { template <class SrcIterator, class SrcAccessor, class DestIteratorX, class DestAccessorX, class DestIteratorY, class DestAccessorY> inline void gaussianGradient(triple<SrcIterator, SrcIterator, SrcAccessor> src, pair<DestIteratorX, DestAccessorX> destx, pair<DestIteratorY, DestAccessorY> desty, double scale); } Usage: #include "vigra/convolution.hxx"
vigra::FImage src(w,h), gradx(w,h), grady(w,h);
...
// calculate gradient vector at scale = 3.0
vigra::gaussianGradient(srcImageRange(src),
destImage(gradx), destImage(grady), 3.0);
|
|
void gaussianSharpening (...) |
|
Perform sharpening function with gaussian filter. This function use the gaussianSmoothing() at first and scale the source image ( src scale tmp dest = (1 + sharpening_factor)*src - sharpening_factor*tmp Preconditions: 1. sharpening_factor >= 0
2. scale >= 0
Declarations: pass arguments explicitly: namespace vigra { template <class SrcIterator, class SrcAccessor, class DestIterator, class DestAccessor> void gaussianSharpening(SrcIterator src_ul, SrcIterator src_lr, SrcAccessor src_acc, DestIterator dest_ul, DestAccessor dest_acc, double sharpening_factor, double scale) } use argument objects in conjuction with Argument Object Factories: namespace vigra { template <class SrcIterator, class SrcAccessor, class DestIterator, class DestAccessor> void gaussianSharpening(triple<SrcIterator, SrcIterator, SrcAccessor> src, pair<DestIterator, DestAccessor> dest, double sharpening_factor, double scale) } Usage: #include "vigra/convolution.hxx"
vigra::FImage src(w,h), dest(w,h);
...
// sharpening with sharpening_factor = 3.0
// smoothing with scale = 0.5
vigra::gaussianSmoothing(srcImageRange(src), destImage(dest), 3.0, 0.5);
|
|
void gaussianSmoothing (...) |
|
Perform isotropic Gaussian convolution.
This function is a shorthand for the concatenation of a call to separableConvolveX() and separableConvolveY() with a Gaussian kernel of the given scale. The function uses Declarations: pass arguments explicitly: namespace vigra { template <class SrcIterator, class SrcAccessor, class DestIterator, class DestAccessor> void gaussianSmoothing(SrcIterator supperleft, SrcIterator slowerright, SrcAccessor sa, DestIterator dupperleft, DestAccessor da, double scale); } use argument objects in conjuction with Argument Object Factories: namespace vigra { template <class SrcIterator, class SrcAccessor, class DestIterator, class DestAccessor> inline void gaussianSmoothing(triple<SrcIterator, SrcIterator, SrcAccessor> src, pair<DestIterator, DestAccessor> dest, double scale); } Usage: #include "vigra/convolution.hxx"
vigra::FImage src(w,h), dest(w,h);
...
// smooth with scale = 3.0
vigra::gaussianSmoothing(srcImageRange(src), destImage(dest), 3.0);
|
|
void hessianMatrixOfGaussian (...) |
|
Filter image with the 2nd derivatives of the Gaussian at the given scale to get the Hessian matrix. The Hessian matrix is a symmetric matrix defined as:
where Declarations: pass arguments explicitly: namespace vigra { template <class SrcIterator, class SrcAccessor, class DestIteratorX, class DestAccessorX, class DestIteratorXY, class DestAccessorXY, class DestIteratorY, class DestAccessorY> void hessianMatrixOfGaussian(SrcIterator supperleft, SrcIterator slowerright, SrcAccessor sa, DestIteratorX dupperleftx, DestAccessorX dax, DestIteratorXY dupperleftxy, DestAccessorXY daxy, DestIteratorY dupperlefty, DestAccessorY day, double scale); } use argument objects in conjuction with Argument Object Factories: namespace vigra { template <class SrcIterator, class SrcAccessor, class DestIteratorX, class DestAccessorX, class DestIteratorXY, class DestAccessorXY, class DestIteratorY, class DestAccessorY> inline void hessianMatrixOfGaussian(triple<SrcIterator, SrcIterator, SrcAccessor> src, pair<DestIteratorX, DestAccessorX> destx, pair<DestIteratorXY, DestAccessorXY> destxy, pair<DestIteratorY, DestAccessorY> desty, double scale); } Usage: #include "vigra/convolution.hxx"
vigra::FImage src(w,h), hxx(w,h), hxy(w,h), hyy(w,h);
...
// calculate Hessian of Gaussian at scale = 3.0
vigra::hessianMatrixOfGaussian(srcImageRange(src),
destImage(hxx), destImage(hxy), destImage(hyy), 3.0);
|
|
void laplacianOfGaussian (...) |
|
Filter image with the Laplacian of Gaussian operator at the given scale. This function calls separableConvolveX() and separableConvolveY() with the appropriate 2nd derivative of Gaussian kernels in x- and y-direction and then sums the results to get the Laplacian. Declarations: pass arguments explicitly: namespace vigra { template <class SrcIterator, class SrcAccessor, class DestIterator, class DestAccessor> void laplacianOfGaussian(SrcIterator supperleft, SrcIterator slowerright, SrcAccessor sa, DestIterator dupperleft, DestAccessor da, double scale); } use argument objects in conjuction with Argument Object Factories: namespace vigra { template <class SrcIterator, class SrcAccessor, class DestIterator, class DestAccessor> inline void laplacianOfGaussian(triple<SrcIterator, SrcIterator, SrcAccessor> src, pair<DestIterator, DestAccessor> dest, double scale); } Usage: #include "vigra/convolution.hxx"
vigra::FImage src(w,h), dest(w,h);
...
// calculate Laplacian of Gaussian at scale = 3.0
vigra::laplacianOfGaussian(srcImageRange(src), destImage(dest), 3.0);
|
|
void simpleSharpening (...) |
|
Perform simple sharpening function. This function use convolveImage( ) with following filter:
-sharpening_factor/16.0, -sharpening_factor/8.0, -sharpening_factor/16.0,
-sharpening_factor/8.0, 1.0+sharpening_factor*0.75, -sharpening_factor/8.0,
-sharpening_factor/16.0, -sharpening_factor/8.0, -sharpening_factor/16.0;
and use Preconditions: 1. sharpening_factor >= 0
2. scale >= 0
Declarations: Declarations: pass arguments explicitly: namespace vigra { template <class SrcIterator, class SrcAccessor, class DestIterator, class DestAccessor> void simpleSharpening(SrcIterator src_ul, SrcIterator src_lr, SrcAccessor src_acc, DestIterator dest_ul, DestAccessor dest_acc, double sharpening_factor) } use argument objects in conjuction with Argument Object Factories: namespace vigra { template <class SrcIterator, class SrcAccessor, class DestIterator, class DestAccessor> inline void simpleSharpening(triple<SrcIterator, SrcIterator, SrcAccessor> src, pair<DestIterator, DestAccessor> dest, double sharpening_factor) { simpleSharpening(src.first, src.second, src.third, dest.first, dest.second, sharpening_factor); } } Usage: #include "vigra/convolution.hxx"
vigra::FImage src(w,h), dest(w,h);
...
// sharpening with sharpening_factor = 0.1
vigra::simpleSharpening(srcImageRange(src), destImage(dest), 0.1);
|
|
void structureTensor (...) |
|
Calculate the Structure Tensor for each pixel of and image, using Gaussian (derivative) filters. The Structure Tensor is is a smoothed version of the Euclidean product of the gradient vector with itself. I.e. it's a symmetric matrix defined as:
where Declarations: pass arguments explicitly: namespace vigra { template <class SrcIterator, class SrcAccessor, class DestIteratorX, class DestAccessorX, class DestIteratorXY, class DestAccessorXY, class DestIteratorY, class DestAccessorY> void structureTensor(SrcIterator supperleft, SrcIterator slowerright, SrcAccessor sa, DestIteratorX dupperleftx, DestAccessorX dax, DestIteratorXY dupperleftxy, DestAccessorXY daxy, DestIteratorY dupperlefty, DestAccessorY day, double inner_scale, double outer_scale); } use argument objects in conjuction with Argument Object Factories: namespace vigra { template <class SrcIterator, class SrcAccessor, class DestIteratorX, class DestAccessorX, class DestIteratorXY, class DestAccessorXY, class DestIteratorY, class DestAccessorY> inline void structureTensor(triple<SrcIterator, SrcIterator, SrcAccessor> src, pair<DestIteratorX, DestAccessorX> destx, pair<DestIteratorXY, DestAccessorXY> destxy, pair<DestIteratorY, DestAccessorY> desty, double nner_scale, double outer_scale); } Usage: #include "vigra/convolution.hxx"
vigra::FImage src(w,h), stxx(w,h), stxy(w,h), styy(w,h);
...
// calculate Structure Tensor at inner scale = 1.0 and outer scale = 3.0
vigra::structureTensor(srcImageRange(src),
destImage(stxx), destImage(stxy), destImage(styy), 1.0, 3.0);
|
|
© Ullrich Köthe (koethe@informatik.uni-hamburg.de) |
html generated using doxygen and Python
|