Create Sobel Kernel#

Synopsis#

Create the Sobel kernel.

Results#

Output:

Neighborhood:
Radius:[1, 1]
Size:[3, 3]
DataBuffer:NeighborhoodAllocator { this = 0x7ffee3e84a00, begin = 0x7fca93256b10, size=9 }
-1
0
1
-2
0
2
-1
0
1

Code#

Python#

#!/usr/bin/env python

import itk

sobelOperator = itk.SobelOperator[itk.F, 2]()
sobelOperator.SetDirection(0)  # Create the operator for the X axis derivative
radius = itk.Size[2]()
radius.Fill(1)
sobelOperator.CreateToRadius(radius)

print(sobelOperator)

for i in range(9):
    print(sobelOperator.GetElement(i))

C++#

#include <itkSobelOperator.h>

int
main()
{
  using SobelOperatorType = itk::SobelOperator<float, 2>;
  SobelOperatorType sobelOperator;
  sobelOperator.SetDirection(0); // Create the operator for the X axis derivative
  itk::Size<2> radius;
  radius.Fill(1);
  sobelOperator.CreateToRadius(radius);

  std::cout << sobelOperator << std::endl;

  for (unsigned int i = 0; i < 9; ++i)
  {
    std::cout << sobelOperator.GetElement(i) << std::endl;
  }
  return EXIT_SUCCESS;
}

Classes demonstrated#

template<typename TPixel, unsigned int VDimension = 2, typename TAllocator = NeighborhoodAllocator<TPixel>>
class SobelOperator : public itk::NeighborhoodOperator<TPixel, VDimension, TAllocator>

A NeighborhoodOperator for performing a directional Sobel edge-detection operation at a pixel location.

SobelOperator is a directional NeighborhoodOperator that should be applied a NeighborhoodIterator using the NeighborhoodInnerProduct method. To create the operator:

1) Set the direction by calling

SetDirection 
2) call
itk::Size<2> radius;
radius.Fill(1);
sobelOperator.CreateToRadius(radius);
3) You may optionally scale the coefficients of this operator using the
ScaleCoefficients 
method. This is useful if you want to take the spacing of the image into account when computing the edge strength. Apply the scaling only after calling to
CreateToRadius 
.

The Sobel Operator in vertical direction for 2 dimensions is

*             -1  -2  -1
*             0    0   0
*             1    2   1
*
*
The Sobel Operator in horizontal direction is for 2 dimensions is
*             -1   0   1
*             -2   0   2
*             -1   0   1
*

The current implementation of the Sobel operator is for 2 and 3 dimensions only. The ND version is planned for future releases.

The extension to 3D is from the publication “Irwin Sobel. An Isotropic 3x3x3 Volume Gradient Operator.

Technical report, Hewlett-Packard Laboratories, April 1995.”

The Sobel operator in 3D has the kernel

* -1 -3 -1   0 0 0  1 3 1
* -3 -6 -3   0 0 0  3 6 3
* -1 -3 -1   0 0 0  1 3 1
*
*    x-1       x     x+1
*

The x kernel is just rotated as required to obtain the kernel in the y and z directions.

Note

SobelOperator does not have any user-declared “special member function”, following the C++ Rule of Zero: the compiler will generate them if necessary.

See

NeighborhoodOperator

See

Neighborhood

See

ForwardDifferenceOperator

See

BackwardDifferenceOperator

ITK Sphinx Examples:

See itk::SobelOperator for additional documentation.