Filtering functions
ImageFiltering.imfilter
— Function.imfilter([T], img, kernel, [border="replicate"], [alg]) --> imgfilt
imfilter([r], img, kernel, [border="replicate"], [alg]) --> imgfilt
imfilter(r, T, img, kernel, [border="replicate"], [alg]) --> imgfilt
Filter an array img
with kernel kernel
by computing their correlation.
kernel[0,0,..]
corresponds to the origin (zero displacement) of the kernel; you can use centered
to place the origin at the array center, or use the OffsetArrays package to set kernel
's indices manually. For example, to filter with a random centered 3x3 kernel, you could use either of the following:
kernel = centered(rand(3,3))
kernel = OffsetArray(rand(3,3), -1:1, -1:1)
kernel
can be specified as an array or as a "factored kernel," a tuple (filt1, filt2, ...)
of filters to apply along each axis of the image. In cases where you know your kernel is separable, this format can speed processing. Each of these should have the same dimensionality as the image itself, and be shaped in a manner that indicates the filtering axis, e.g., a 3x1 filter for filtering the first dimension and a 1x3 filter for filtering the second dimension. In two dimensions, any kernel passed as a single matrix is checked for separability; if you want to eliminate that check, pass the kernel as a single-element tuple, (kernel,)
.
Optionally specify the border
, as one of Fill(value)
, "replicate"
, "circular"
, "symmetric"
, "reflect"
, NA()
, or Inner()
. The default is "replicate"
. These choices specify the boundary conditions, and therefore affect the result at the edges of the image. See padarray
for more information.
alg
allows you to choose the particular algorithm: FIR()
(finite impulse response, aka traditional digital filtering) or FFT()
(Fourier-based filtering). If no choice is specified, one will be chosen based on the size of the image and kernel in a way that strives to deliver good performance. Alternatively you can use a custom filter type, like KernelFactors.IIRGaussian
.
Optionally, you can control the element type of the output image by passing in a type T
as the first argument.
You can also dispatch to different implementations by passing in a resource r
as defined by the ComputationalResources package. For example,
imfilter(ArrayFire(), img, kernel)
would request that the computation be performed on the GPU using the ArrayFire libraries.
See also: imfilter!
, centered
, padarray
, Pad
, Fill
, Inner
, KernelFactors.IIRGaussian
.
ImageFiltering.imfilter!
— Function.imfilter!(imgfilt, img, kernel, [border="replicate"], [alg])
imfilter!(r, imgfilt, img, kernel, border, [inds])
imfilter!(r, imgfilt, img, kernel, border::NoPad, [inds=indices(imgfilt)])
Filter an array img
with kernel kernel
by computing their correlation, storing the result in imgfilt
.
The indices of imgfilt
determine the region over which the filtered image is computed–-you can use this fact to select just a specific region of interest, although be aware that the input img
might still get padded. Alteratively, explicitly provide the indices inds
of imgfilt
that you want to calculate, and use NoPad
boundary conditions. In such cases, you are responsible for supplying appropriate padding: img
must be indexable for all of the locations needed for calculating the output. This syntax is best-supported for FIR filtering; in particular, that that IIR filtering can lead to results that are inconsistent with respect to filtering the entire array.
See also: imfilter
.
ImageFiltering.imgradients
— Function.imgradients(img, kernelfun=KernelFactors.ando3, border="replicate") -> gimg1, gimg2, ...
Estimate the gradient of img
at all points of the image, using a kernel specified by kernelfun
. The gradient is returned as a tuple-of-arrays, one for each dimension of the input; gimg1
corresponds to the derivative with respect to the first dimension, gimg2
to the second, and so on. At the image edges, border
is used to specify the boundary conditions.
kernelfun
may be one of the filters defined in the KernelFactors
module, or more generally any function which satisfies the following interface:
kernelfun(extended::NTuple{N,Bool}, d) -> kern_d
kern_d
is the kernel for producing the derivative with respect to the d
th dimension of an N
-dimensional array. extended[i]
is true if the image is of size > 1 along dimension i
. kern_d
may be provided as a dense or factored kernel, with factored representations recommended when the kernel is separable.
imgradients(img, points, [kernelfunc], [border]) -> G
Performs edge detection filtering in the N-dimensional array img
. Gradients are computed at specified points
(or indexes) in the array.
All kernel functions are specified as KernelFactors.func
. For 2d images, the choices for func
include sobel
, prewitt
, ando3
, ando4
, and ando5
. For other dimensionalities, the ando4
and ando5
kernels are not available.
Border options:"replicate"
, "circular"
, "reflect"
, "symmetric"
.
Returns a 2D array G
with the gradients as rows. The number of rows is the number of points at which the gradient was computed and the number of columns is the dimensionality of the array.
ImageFiltering.MapWindow.mapwindow
— Function.mapwindow(f, img, window, [border="replicate"]) -> imgf
Apply f
to sliding windows of img
, with window size or indices specified by window
. For example, mapwindow(median!, img, window)
returns an Array
of values similar to img
(median-filtered, of course), whereas mapwindow(extrema, img, window)
returns an Array
of (min,max)
tuples over a window of size window
centered on each point of img
.
The function f
receives a buffer buf
for the window of data surrounding the current point. If window
is specified as a Dims-tuple (tuple-of-integers), then all the integers must be odd and the window is centered around the current image point. For example, if window=(3,3)
, then f
will receive an Array buf
corresponding to offsets (-1:1, -1:1)
from the imgf[i,j]
for which this is currently being computed. Alternatively, window
can be a tuple of AbstractUnitRanges, in which case the specified ranges are used for buf
; this allows you to use asymmetric windows if needed.
border
specifies how the edges of img
should be handled; see imfilter
for details.
For functions that can only take AbstractVector
inputs, you might have to first specialize default_shape
:
f = v->quantile(v, 0.75)
ImageFiltering.MapWindow.default_shape(::typeof(f)) = vec
and then mapwindow(f, img, (m,n))
should filter at the 75th quantile.
See also: imfilter
.
Kernel
ImageFiltering.Kernel.sobel
— Function.diff1, diff2 = sobel()
Return kernels for two-dimensional gradient compution using the Sobel operator. diff1
computes the gradient along the first (y) dimension, and diff2
computes the gradient along the second (x) dimension.
See also: KernelFactors.sobel
, Kernel.prewitt
, Kernel.ando3
.
ImageFiltering.Kernel.prewitt
— Function.diff1, diff2 = prewitt()
Return kernels for two-dimensional gradient compution using the Prewitt operator. diff1
computes the gradient along the first (y) dimension, and diff2
computes the gradient along the second (x) dimension.
See also: KernelFactors.prewitt
, Kernel.sobel
, Kernel.ando3
.
ImageFiltering.Kernel.ando3
— Function.diff1, diff2 = ando3()
Return 3x3 kernels for two-dimensional gradient compution using the optimal "Ando" filters. diff1
computes the gradient along the y-axis (first dimension), and diff2
computes the gradient along the x-axis (second dimension).
Citation
Ando Shigeru, IEEE Trans. Pat. Anal. Mach. Int., vol. 22 no 3, March 2000
See also: KernelFactors.ando3
, Kernel.ando4
, Kernel.ando5
.
ImageFiltering.Kernel.ando4
— Function.diff1, diff2 = ando4()
Return 4x4 kernels for two-dimensional gradient compution using the optimal "Ando" filters. diff1
computes the gradient along the y-axis (first dimension), and diff2
computes the gradient along the x-axis (second dimension).
Citation
Ando Shigeru, IEEE Trans. Pat. Anal. Mach. Int., vol. 22 no 3, March 2000
See also: KernelFactors.ando4
, Kernel.ando3
, Kernel.ando5
.
ImageFiltering.Kernel.ando5
— Function.diff1, diff2 = ando5()
Return 5x5 kernels for two-dimensional gradient compution using the optimal "Ando" filters. diff1
computes the gradient along the y-axis (first dimension), and diff2
computes the gradient along the x-axis (second dimension).
Citation
Ando Shigeru, IEEE Trans. Pat. Anal. Mach. Int., vol. 22 no 3, March 2000
See also: KernelFactors.ando5
, Kernel.ando3
, Kernel.ando4
.
ImageFiltering.Kernel.gaussian
— Function.gaussian((σ1, σ2, ...), [(l1, l2, ...]) -> g
gaussian(σ) -> g
Construct a multidimensional gaussian filter, with standard deviation σd
along dimension d
. Optionally provide the kernel length l
, which must be a tuple of the same length.
If σ
is supplied as a single number, a symmetric 2d kernel is constructed.
See also: KernelFactors.gaussian
.
ImageFiltering.Kernel.DoG
— Function.DoG((σp1, σp2, ...), (σm1, σm2, ...), [l1, l2, ...]) -> k
DoG((σ1, σ2, ...)) -> k
DoG(σ::Real) -> k
Construct a multidimensional difference-of-gaussian kernel k
, equal to gaussian(σp, l)-gaussian(σm, l)
. When only a single σ
is supplied, the default is to choose σp = σ, σm = √2 σ
. Optionally provide the kernel length l
; the default is to extend by two max(σp,σm)
in each direction from the center. l
must be odd.
If σ
is provided as a single number, a symmetric 2d DoG kernel is returned.
See also: KernelFactors.IIRGaussian
.
ImageFiltering.Kernel.LoG
— Function.LoG((σ1, σ2, ...)) -> k
LoG(σ) -> k
Construct a Laplacian-of-Gaussian kernel k
. σd
is the gaussian width along dimension d
. If σ
is supplied as a single number, a symmetric 2d kernel is returned.
See also: KernelFactors.IIRGaussian
and Kernel.Laplacian
.
ImageFiltering.Kernel.Laplacian
— Type.Laplacian((true,true,false,...))
Laplacian(dims, N)
Lacplacian()
Laplacian kernel in N
dimensions, taking derivatives along the directions marked as true
in the supplied tuple. Alternatively, one can pass dims
, a listing of the dimensions for differentiation. (However, this variant is not inferrable.)
Laplacian()
is the 2d laplacian, equivalent to Laplacian((true,true))
.
The kernel is represented as an opaque type, but you can use convert(AbstractArray, L)
to convert it into array format.
KernelFactors
ImageFiltering.KernelFactors.sobel
— Function.kern1, kern2 = sobel()
Factored Sobel filters for dimensions 1 and 2 of a two-dimensional image. Each is a 2-tuple of one-dimensional filters.
kern = sobel(extended::NTuple{N,Bool}, d)
Return a factored Sobel filter for computing the gradient in N
dimensions along axis d
. If extended[dim]
is false, kern
will have size 1 along that dimension.
ImageFiltering.KernelFactors.prewitt
— Function.kern1, kern2 = prewitt()
returns factored Prewitt filters for dimensions 1 and 2 of your image
kern = prewitt(extended::NTuple{N,Bool}, d)
Return a factored Prewitt filter for computing the gradient in N
dimensions along axis d
. If extended[dim]
is false, kern
will have size 1 along that dimension.
ImageFiltering.KernelFactors.ando3
— Function.kern1, kern2 = ando3()
returns optimal 3x3 gradient filters for dimensions 1 and 2 of your image, as defined in Ando Shigeru, IEEE Trans. Pat. Anal. Mach. Int., vol. 22 no 3, March 2000.
See also: Kernel.ando3
, KernelFactors.ando4
, KernelFactors.ando5
.
kern = ando3(extended::NTuple{N,Bool}, d)
Return a factored Ando filter (size 3) for computing the gradient in N
dimensions along axis d
. If extended[dim]
is false, kern
will have size 1 along that dimension.
ImageFiltering.KernelFactors.ando4
— Function.kern1, kern2 = ando4()
returns separable approximations of the optimal 4x4 filters for dimensions 1 and 2 of your image, as defined in Ando Shigeru, IEEE Trans. Pat. Anal. Mach. Int., vol. 22 no 3, March 2000.
See also: Kernel.ando4
.
kern = ando4(extended::NTuple{N,Bool}, d)
Return a factored Ando filter (size 4) for computing the gradient in N
dimensions along axis d
. If extended[dim]
is false, kern
will have size 1 along that dimension.
ImageFiltering.KernelFactors.ando5
— Function.kern1, kern2 = ando5_sep()
returns separable approximations of the optimal 5x5 gradient filters for dimensions 1 and 2 of your image, as defined in Ando Shigeru, IEEE Trans. Pat. Anal. Mach. Int., vol. 22 no 3, March 2000.
See also: Kernel.ando5
.
kern = ando5(extended::NTuple{N,Bool}, d)
Return a factored Ando filter (size 5) for computing the gradient in N
dimensions along axis d
. If extended[dim]
is false, kern
will have size 1 along that dimension.
ImageFiltering.KernelFactors.gaussian
— Function.gaussian(σ::Real, [l]) -> g
Construct a 1d gaussian kernel g
with standard deviation σ
, optionally providing the kernel length l
. The default is to extend by two σ
in each direction from the center. l
must be odd.
gaussian((σ1, σ2, ...), [l]) -> (g1, g2, ...)
Construct a multidimensional gaussian filter as a product of single-dimension factors, with standard deviation σd
along dimension d
. Optionally provide the kernel length l
, which must be a tuple of the same length.
ImageFiltering.KernelFactors.IIRGaussian
— Function.IIRGaussian([T], σ; emit_warning::Bool=true)
Construct an infinite impulse response (IIR) approximation to a Gaussian of standard deviation σ
. σ
may either be a single real number or a tuple of numbers; in the latter case, a tuple of such filters will be created, each for filtering a different dimension of an array.
Optionally specify the type T
for the filter coefficients; if not supplied, it will match σ
(unless σ
is not floating-point, in which case Float64
will be chosen).
Citation
I. T. Young, L. J. van Vliet, and M. van Ginkel, "Recursive Gabor Filtering". IEEE Trans. Sig. Proc., 50: 2798-2805 (2002).
TriggsSdika(a, b, scale, M)
Defines a kernel for one-dimensional infinite impulse response (IIR) filtering. a
is a "forward" filter, b
a "backward" filter, M
is a matrix for matching boundary conditions at the right edge, and scale
is a constant scaling applied to each element at the conclusion of filtering.
Citation
B. Triggs and M. Sdika, "Boundary conditions for Young-van Vliet recursive filtering". IEEE Trans. on Sig. Proc. 54: 2365-2367 (2006).
TriggsSdika(ab, scale)
Create a symmetric Triggs-Sdika filter (with a = b = ab
). M
is calculated for you. Only length 3 filters are currently supported.
Kernel utilities
ImageFiltering.centered
— Function.centered(kernel) -> shiftedkernel
Shift the origin-of-coordinates to the center of kernel
. The center-element of kernel
will be accessed by shiftedkernel[0, 0, ...]
.
This function makes it easy to supply kernels using regular Arrays, and provides compatibility with other languages that do not support arbitrary indices.
See also: imfilter
.
ImageFiltering.KernelFactors.kernelfactors
— Function.kernelfactors(factors::Tuple)
Prepare a factored kernel for filtering. If passed a 2-tuple of vectors of lengths m
and n
, this will return a 2-tuple of ReshapedVector
s that are effectively of sizes m×1
and 1×n
. In general, each successive factor
will be reshaped to extend along the corresponding dimension.
If passed a tuple of general arrays, it is assumed that each is shaped appropriately along its "leading" dimensions; the dimensionality of each is "extended" to N = length(factors)
, appending 1s to the size as needed.
ImageFiltering.Kernel.reflect
— Function.reflect(kernel) --> reflectedkernel
Compute the pointwise reflection around 0, 0, ... of the kernel kernel
. Using imfilter
with a reflectedkernel
performs convolution, rather than correlation, with respect to the original kernel
.
Boundaries and padding
ImageFiltering.padarray
— Function.padarray([T], img, border) --> imgpadded
Generate a padded image from an array img
and a specification border
of the boundary conditions and amount of padding to add. border
can be a Pad
, Fill
, or Inner
object.
Optionally provide the element type T
of imgpadded
.
ImageFiltering.Pad
— Type.Pad
is a type that stores choices about padding. Instances must set style
, a Symbol specifying the boundary conditions of the image, one of:
:replicate
(repeat edge values to infinity):circular
(image edges "wrap around"):symmetric
(the image reflects relative to a position between pixels):reflect
(the image reflects relative to the edge itself)
The default value is :replicate
.
It's worth emphasizing that padding is most straightforwardly specified as a string,
imfilter(img, kernel, "replicate")
rather than
imfilter(img, kernel, Pad(:replicate))
ImageFiltering.Fill
— Type.Fill(val)
Fill(val, lo, hi)
Pad the edges of the image with a constant value, val
.
Optionally supply the extent of the padding, see Pad
.
Example:
imfilter(img, kernel, Fill(zero(eltype(img))))
ImageFiltering.Inner
— Type.Inner()
Inner(lo, hi)
Indicate that edges are to be discarded in filtering, only the interior of the result it to be returned.
Example:
imfilter(img, kernel, Inner())
ImageFiltering.NA
— Type.NA()
NA(lo, hi)
Choose filtering using "NA" (Not Available) boundary conditions. This is most appropriate for filters that have only positive weights, such as blurring filters. Effectively, the output pixel value is normalized in the following way:
filtered img with Fill(0) boundary conditions
output = ---------------------------------------------
filtered 1 with Fill(0) boundary conditions
As a consequence, filtering has the same behavior as nanmean
. Indeed, invalid pixels in img
can be marked as NaN
and then they are effectively omitted from the filtered result.
ImageFiltering.NoPad
— Type.NoPad()
NoPad(border)
Indicates that no padding should be applied to the input array, or that you have already pre-padded the input image. Passing a border
object allows you to preserve "memory" of a border choice; it can be retrieved by indexing with []
.
Example
np = NoPad(Pad(:replicate))
imfilter!(out, img, kernel, np)
runs filtering directly, skipping any padding steps. Every entry of out
must be computable using in-bounds operations on img
and kernel
.
Algorithms
ImageFiltering.Algorithm.FIR
— Type.Filter using a direct algorithm
ImageFiltering.Algorithm.FFT
— Type.Filter using the Fast Fourier Transform
ImageFiltering.Algorithm.IIR
— Type.Filter with an Infinite Impulse Response filter
ImageFiltering.Algorithm.Mixed
— Type.Filter with a cascade of mixed types (IIR, FIR)
Internal machinery
ReshapedOneD{N,Npre}(data)
Return an object of dimensionality N
, where data
must have dimensionality 1. The indices are 0:0
for the first Npre
dimensions, have the indices of data
for dimension Npre+1
, and are 0:0
for the remaining dimensions.
data
must support eltype
and ndims
, but does not have to be an AbstractArray.
ReshapedOneDs allow one to specify a "filtering dimension" for a 1-dimensional filter.