Reference

The following summarizes the documentation of types and methods provided by the LocalFilters package. This information is also available from the REPL by typing ? followed by the name of a method or a type.

Index

Simple linear filters

LocalFilters provides a number of shift-invariant linear filters.

LocalFilters.correlateFunction
correlate(A, B=3) -> dst

yields the discrete correlation of the array A by the kernel defined by B. The result dst is an array similar to A.

Using Idx(A) to denote the set of valid indices for array A and assuming B is an array of numerical values, the discrete convolution of A by B writes:

dst = similar(A, T)
for i ∈ Idx(dst)
    v = zero(T)
    @inbounds for j ∈ Idx(A) ∩ (i + Idx(B))
        v += A[j]*B[j-i]
    end
    dst[i] = v
end

with T the type of the sum of the products of the elements of A and B, and where Idx(A) ∩ (i + Idx(B)) denotes the subset of indices j such that j ∈ Idx(A) and j - i ∈ Idx(B) and thus for which A[j] and B[j-i] are valid.

See also correlate! and convolve.

source
LocalFilters.convolveFunction
convolve(A, B=3)

yields the discrete convolution of array A by the kernel defined by B. The result dst is an array similar to A.

Using Idx(A) to denote the set of valid indices for array A and assuming B is an array of values, the discrete convolution of A by B writes:

dst = similar(A, T)
for i ∈ Idx(dst)
    v = zero(T)
    @inbounds for j ∈ Idx(A) ∩ (i - Idx(B))
        v += A[j]*B[i-j]
    end
    dst[i] = v
end

with T the type of the sum of the products of the elements of A and B, and where Idx(A) ∩ (i - Idx(B)) denotes the subset of indices j such that j ∈ Idx(A) and i - j ∈ Idx(B) and thus for which A[j] and B[i-j] are valid.

See also convolve! and localfilter!.

source
LocalFilters.localmeanFunction
localmean([T,] A, B=3; null=zero(T), order=FORWARD_FILTER)

yields the local mean of A in a neighborhood defined by B. The result is an array similar to A. If B is not specified, the neighborhood is a hyper-rectangular sliding window of size 3 in every dimension. Otherwise, B may be specified as a Cartesian box, or as an array of booleans of same number of dimensions as A. If B is a single odd integer (as it is by default), the neighborhood is assumed to be a hyper-rectangular sliding window of size B in every dimension.

Optional argument T is to specify the element type of the result.

Keyword null may be used to specify the value of the result where the sum of the weights in a local neighborhood is zero. By default, null = zero(T).

Keyword order specifies the filter direction, FORWARD_FILTER by default.

See also localmean! and localfilter!.

source
LocalFilters.localmean!Function
localmean!(dst, A, B=3; null=zero(eltype(dst)), order=FORWARD_FILTER) -> dst

overwrites dst with the local mean of A in a neighborhood defined by B and returns dst.

Keyword null may be used to specify the value of the result where the sum of the weights in the a neighborhood is zero.

Keyword order specifies the filter direction, FORWARD_FILTER by default.

See also localmean and localfilter!.

source

Mathematical morphology

LocalFilters.erodeFunction
erode(A, B=3; order=FORWARD_FILTER, slow=false) -> Amin

yields the erosion of A by the structuring element defined by B. The returned result, Amin, is similar to A (same size and type) and its values are the local minima of A in the neighborhood defined by B.

Keyword order specifies the filter direction, FORWARD_FILTER by default.

If B is not an N-dimensional array, kernel(Dims{N},B) is called to build a kernel with N = ndims(A) the number of dimensions of A.

If the structuring element B is equivalent to a simple hyper-rectangular sliding window (which is the case by default) and unless keyword slow is true, the much faster van Herk-Gil-Werman algorithm is used.

An erosion is one of the most basic operations of mathematical morphology. See erode! for an in-place version of the method, dilate for retrieving the local maxima, and localextrema for performing an erosion and a dilation in a single pass.

source
LocalFilters.erode!Function
erode!(Amin, A, B=3; order=FORWARD_FILTER, slow=false) -> Amin

overwrites Amin with the erosion of the array A by the structuring element defined by B and returns Amin.

Keyword order specifies the filter direction, FORWARD_FILTER by default.

If the structuring element B is equivalent to a simple hyper-rectangular sliding window (which is the case by default) and unless keyword slow is true, the much faster van Herk-Gil-Werman algorithm is used and the operation can be done in-place. That is, A and Amin can be the same arrays. In that case, the following syntax is allowed:

erode!(A, B=3; order=FORWARD_FILTER, ) -> A

See erode for an out-of-place version and for more information.

source
LocalFilters.dilateFunction
dilate(A, B=3; order=FORWARD_FILTER, slow=false) -> Amax

yields the dilation of A by the structuring element defined by B. The returned result, Amax, is similar to A (same size and type) and its values are the local maxima of A in the neighborhood defined by B.

Keyword order specifies the filter direction, FORWARD_FILTER by default.

If B is not an N-dimensional array, kernel(Dims{N},B) is called to build a kernel with N = ndims(A) the number of dimensions of A.

If the structuring element B is equivalent to a simple hyper-rectangular sliding window (which is the case by default) and unless keyword slow is true, the much faster van Herk-Gil-Werman algorithm is used.

A dilation is one of the most basic operations of mathematical morphology. See dilate! for an in-place version of the method, erode for retrieving the local minima, and localextrema for performing an erosion and a dilation in a single pass.

source
LocalFilters.dilate!Function
dilate!(Amax, A, B=3; order=FORWARD_FILTER, slow=false) -> Amax

overwrites Amax with a dilation of the array A by the structuring element defined by B and returns Amax.

Keyword order specifies the filter direction, FORWARD_FILTER by default.

If the structuring element B is equivalent to a simple hyper-rectangular sliding window (which is the case by default) and unless keyword slow is true, the much faster van Herk-Gil-Werman algorithm is used and the operation can be done in-place. That is, A and Amin can be the same arrays. In that case, the following syntax is allowed:

dilate!(A, B=3; order=FORWARD_FILTER) -> A

See dilate for an out-of-place version and for more information.

source
LocalFilters.localextremaFunction
localextrema(A, B=3; order=FORWARD_FILTER) -> Amin, Amax

yields the results of performing an erosion and a dilation of A by the structuring element defined by B in a single pass. Calling this method is usually almost twice as fast as calling erode and dilate.

Keyword order specifies the filter direction, FORWARD_FILTER by default.

See localextrema! for an in-place version of the method, and erode or dilate for a description of these operations.

source
LocalFilters.localextrema!Function
localextrema!(Amin, Amax, A, B=3; order=FORWARD_FILTER) -> Amin, Amax

overwrites Amin and Amax with, respectively, an erosion and a dilation of the array A by the structuring element defined by B in a single pass.

Keyword order specifies the filter direction, FORWARD_FILTER by default.

See localextrema for an out-of-place version for more information.

source
LocalFilters.closingFunction
closing(A, B=3; order=FORWARD_FILTER, slow=false) -> dst

yields a closing of array A by the structuring element defined by B. A closing is a dilation followed by an erosion. The result dst is an array similar to A.

See closing! for an in-place version of the method, opening for a related filter, and erode or dilate for a description of these operations and for the meaning of keywords.

source
LocalFilters.closing!Function
closing!(dst, wrk, A, B=3; order=FORWARD_FILTER, slow=false) -> dst

overwrites dst with the result of a closing of A by the structuring element defined by B using wrk as a workspace array. The arguments dst, wrk, and A must be similar arrays, dst and A may be identical, but wrk must not be the same array as A or dst. The destination dst is returned.

See closing for a description of this kind of filter and for the meaning of the arguments and keywords.

source
LocalFilters.openingFunction
opening(A, B=3; order=FORWARD_FILTER, slow=false) -> dst

yields an opening of array A by the structuring element defined by B. An opening is an erosion followed by a dilation. The result dst is an array similar to A.

See opening! for an in-place version of the method, closing for a related filter, and erode or dilate for a description of these operations and for the meaning of keywords.

source
LocalFilters.opening!Function
opening!(dst, wrk, A, B=3; order=FORWARD_FILTER, slow=false) -> dst

overwrites dst with the result of an opening of A by the structuring element defined by B using wrk as a workspace array. The arguments dst, wrk, and A must be similar arrays, dst and A may be identical, but wrk must not be the same array as A or dst. The destination dst is returned.

See opening for a description of this kind of filter and for the meaning of the arguments and keywords.

source
LocalFilters.top_hatFunction
top_hat(A, B=3[, C]; order=FORWARD_FILTER, slow=false) -> dst

performs a summit detection by applying a top-hat filter to array A using the structuring element defined by B for the feature detection. Top-hat filtering is equivalent to:

dst = A .- opening(A, B)

Optional argument C specifies the structuring element for smoothing A prior to top-hat filtering. If B and C are specified as the radii of the structuring elements, then C should be smaller than B. For instance:

top_hat(bitmap, 3, 1)

may be used to detect text or lines in a bitmap image.

Keyword order specifies the filter(s) direction(s), FORWARD_FILTER by default. If C is specified, order may be a 2-tuple to specify a first order for B and a second one for C.

See bottom_hat for a related operation, LocalFilters.top_hat! for an in-place version.

source
LocalFilters.top_hat!Function
LocalFilters.top_hat!(dst, wrk, A, B=3; order=FORWARD_FILTER, slow=false) -> dst

overwrites dst with the result of a top-hat filter applied to A with structuring element B, and using wrk as a workspace whose contents is not preserved. The arguments A, dst, and wrk must be similar but different arrays. The destination dst is returned.

See also top_hat for more details.

source
LocalFilters.bottom_hatFunction
bottom_hat(A, B=3[, C]; order=FORWARD_FILTER, slow=false) -> dst

performs a valley detection by applying a bottom-hat filter to array A using the structuring element defined by B for the feature detection. Bottom-hat filtering is equivalent to:

dst = closing(A, B) .- A

Optional argument C specifies the structuring element for smoothing A prior to bottom-hat filtering. If B and C are specified as the radii of the structuring elements, then C should be smaller than B.

Keyword order specifies the filter(s) direction(s), FORWARD_FILTER by default. If C is specified, order may be a 2-tuple to specify a first order for B and a second one for C.

See top_hat for a related operation, LocalFilters.bottom_hat! for an in-place version.

source
LocalFilters.bottom_hat!Function
LocalFilters.bottom_hat!(dst, wrk, A, B=3; order=FORWARD_FILTER, slow=false) -> dst

overwrites dst with the result of a bottom-hat filter applied to A with structuring element B and optional smoothing element C. Argument wrk is a workspace array whose contents is not preserved. The arguments A, dst, and wrk must be similar but different arrays. The destination dst is returned.

See also bottom_hat for more details.

source

Other non-linear filters

LocalFilters.bilateralfilterFunction
bilateralfilter([T = float(eltype(A)),] A, F, G...=3; order = FORWARD_FILTER)

yields the result of applying the bilateral filter on array A.

Argument F specifies how to smooth the differences in values. It can be:

  • a function, say f, which is called as f(A[i],A[j]) to yield a nonnegative weight for i the central index and j the index in a nearby position;

  • a positive real, say σ, which is assumed to be the standard deviation of a Gaussian.

Arguments G, ... specify the settings of the distance filter for smoothing differences in coordinates. There are several possibilities:

  • G... = wgt an array of nonnegative weights or of Booleans. The axes of wgt must have offsets so that the zero index is part of the indices of wgt.

  • G... = f, w with f a function and w any kind of argument that can be used to build a window win specifying the extension of the neighborhood. The value of the distance filter will be max(f(i),0) for all Cartesian index i of win such that win[i] is true. See kernel for the different ways to specify a window.

  • G... = σ or , G... = σ, w with σ a positive real assumed to be the standard deviation of a Gaussian function and w any kind of argument that can be used to build a window win specifying the extension of the neighborhood. If w is not specified, a default window of size ±3σ is assumed.

Optional argument T is to specify the element type of the result. This is needed if the default is unsuitable.

See bilateralfilter! for an in-place version of this function and see Wikipedia for a description of the bilateral filter.

source
LocalFilters.bilateralfilter!Function
bilateralfilter!(dst, A, F, G...; order = FORWARD_FILTER) -> dst

overwrites dst with the result of applying the bilateral filter on array A and returns dst.

See bilateralfilter for a description of the other arguments than dst and see Wikipedia for a description of the bilateral filter.

source

Methods to build local filters

LocalFilters.localfilterFunction
localfilter([T=eltype(A),] A, B, initial, update, final=identity; kwds...) -> dst

out of place version of localfilter! which is equivalent to:

localfilter!(similar(A, T), A, B, initial, update, final=identity; kwds...)

Optional argument T is to specify the element type of the result; by default, T is the element type of A.

source
LocalFilters.localfilter!Function
localfilter!(dst, A, B, initial, update::Function, final::Function=identity;
             order=FORWARD_FILTER) -> dst

overwrites the destination dst with the result of a local filter applied to the source A, on a relative neighborhood defined by B, and implemented by initial, update, and final. The initial argument may be a function or not. The purpose of these latter arguments is explained by the following pseudo-codes implementing the local filtering. Keyword order specifies the filter direction, FORWARD_FILTER by default.

If order = FORWARD_FILTER:

@inbounds for i ∈ indices(dst)
    v = initial isa Function ? initial(A[i]) : initial
    for j ∈ indices(A) ∩ (indices(B) + i)
        v = update(v, A[j], B[j-i])
    end
    dst[i] = final(v)
end

else if order = REVERSE_FILTER:

@inbounds for i ∈ indices(dst)
    v = initial isa Function ? initial(A[i]) : initial
    for j ∈ indices(A) ∩ (i - indices(B))
        v = update(v, A[j], B[i-j])
    end
    dst[i] = final(v)
end

where indices(A) denotes the range of indices of any array A while indices(B) + i and i - indices(B) respectively denote the set of indices j such that j - i ∈ indices(B) and i - j ∈ indices(B). In other words, j ∈ indices(A) ∩ (i - indices(B)) means all indices j such that j ∈ indices(A) and i - j ∈ indices(B) so that A[j] and B[i-j] are in-bounds.

If initial is a function, the initial value of the state variable v in the above pseudo-codes is given by v = initial(A[i]) with i the current index in dst. Hence, in that case, the destination array dst and the source array src must have the same axes.

For example, implementing a local minimum filter (that is, an erosion), is as simple as:

localfilter!(dst, A, B, typemax(eltype(A)),
             (v,a,b) -> ifelse(b, min(v,a), v))

As another example, implementing a convolution by B writes:

T = promote_type(eltype(A), eltype(B))
localfilter!(dst, A, B, zero(T), (v,a,b) -> v + a*b)
source
LocalFilters.localmapFunction
localmap(f, [T=eltype(A),] A, B=3; null=zero(T), order=FORWARD_FILTER)

for each position in A, applies the function f to the values of A extracted from the neighborhood defined by B.

Optional argument T is to specify the element type of the result; by default, T is the element type of A.

Keyword order specifies the filter direction, FORWARD_FILTER by default.

Remarks

  • The function f is never called with an empty vector of values. Keyword null may be used to specify the value of the result where the neighborhood is empty. By default, null = zero(T) with T the element type of the result.

  • The vector of values passed to f may be modified by f if needed (for example for faster sorting of the values).

Examples

With argument f set to minimum or maximum, localmap respectively yields the erosion and the dilation of the input array. However erode and dilate methods are faster until localmap is specialized for these functions.

Applying a median filter of the 2-dimensional image img in a sliding 5×5 window can be done by:

using Statistics
med = localmap(median!, img, 5; eltype=float(eltype(img)), null=NaN)
source
LocalFilters.localmap!Function
localmap!(f, dst, A, B=3; null=zero(eltype(dst)), order=FORWARD_FILTER)

set each entry of dst, to the result of applying the function f to the values of A extracted from the neighborhood defined by B.

The function f is never called with an empty vector of values. Keyword null may be used to specify the value of the result where the neighborhood is empty. By default, null = zero(T) with T the element type of the result.

Keyword order specifies the filter direction, FORWARD_FILTER by default.

source
LocalFilters.localreduceFunction
localreduce(op, [T=eltype(A),] A, dims, rngs; kwds...) -> dst

yields the local reduction by the associative binary operator op of the values of A into contiguous hyper-rectangular neighborhoods defined by the interval(s) rngs along dimension(s) dims of A. The algorithm of van Herk-Gil-Werman is used to compute the reduction.

Optional argument T is to specify the element type of the result.

Argument dims specifies along which dimension(s) of A the local reduction is to be applied, it can be a single integer, a tuple of integers, or a colon : to apply the operation to all dimensions. Dimensions are processed in the order given by dims (the same dimension may appear several times) and there must be a matching interval in rngs to specify the structuring element (except that, if rngs is a single interval, it is used for every dimension in dims). An interval is either an integer or an integer valued unit range in the form kmin:kmax. An interval specified as a single integer is treated as an approximately centered range of this length.

Keyword order specifies the filter direction, FORWARD_FILTER by default.

Keyword work may be used to provide a workspace array of type Vector{T} which is automatically resized as needed.

Assuming a mono-dimensional array A, the single reduction pass:

dst = localreduce(op, A, :, rng)

amounts to computing (assuming forward ordering):

dst[j] =  A[i+kmin] ⋄ A[i+kmin+1] ⋄ ... ⋄ A[i+kmax-1] ⋄ A[i+kmax]

for all j ∈ axes(dst,1), with x ⋄ y = op(x, y), kmin = first(rng) and kmax = last(rng). Note that if kmin = kmax = k, the result of the filter is to operate a simple shift by k along the corresponding dimension and has no effects if k = 0. This can be exploited to not filter some dimension(s).

Flat boundary conditions are assumed for A[i+k] in the above formula.

Examples

The morphological erosion (local minimum) of the array A on a centered structuring element of width 7 in every dimension can be obtained by:

localreduce(min, A, :, -3:3)

or equivalently by:

localreduce(min, A, :, 7)

Index interval 0:0 may be specified to do nothing along the corresponding dimension. For instance, assuming A is a three-dimensional array:

localreduce(max, A, :, (-3:3, 0:0, -4:4))

yields the morphological dilation (i.e. local maximum) of A in a centered local neighborhood of size 7×1×9 (nothing is done along the second dimension). The same result may be obtained with:

localreduce(max, A, (1,3), (-3:3, -4:4))

where the second dimension is omitted from the list of dimensions.

The local average of the two-dimensional array A on a centered sliding window of size 11×11 can be computed as:

localreduce(+, A, :, (-5:5, -5:5)) ./ 11^2

See localreduce! for an in-place version of the method.

source
LocalFilters.localreduce!Function
localreduce!(op, [dst = A,] A, dims, rngs; kwds...)

overwrites the contents of dst with the local reduction by the associative binary operator op of the values of A into contiguous hyper-rectangular neighborhoods defined by the interval(s) rngs along dimension(s) dims of A. Except if a single dimension of interest is specified by dims, the destination dst must have the same indices as the source A (that is, axes(dst) == axes(A)). Operation may be done in-place and dst and A can be the same; this is the default behavior if dst is not specified. The algorithm of van Herk-Gil-Werman is used to compute the reduction.

See localreduce for a full description of the method and for accepted keywords.

The in-place morphological erosion (local minimum) of the array A on a centered structuring element of width 7 in every dimension can be obtained by:

localreduce!(min, A, :, -3:3)
source

Constants

LocalFilters.FORWARD_FILTERConstant
FORWARD_FILTER

is an exported constant object used to indicate forward ordering of indices in local filter operations. It can be called as:

FORWARD_FILTER(i, j) -> j - i

to yield the index in the filter kernel. See also REVERSE_FILTER for reverse ordering and LocalFilters.localindices for building a range of valid indices j.

source
LocalFilters.REVERSE_FILTERConstant
REVERSE_FILTER

is an exported constant object used to indicate reverse ordering of indices in local filter operations. It can be called as:

REVERSE_FILTER(i, j) -> i - j

to yield the index in the filter kernel. See also FORWARD_FILTER for forward ordering and LocalFilters.localindices for building a range of valid indices j.

source

Neighborhoods, boxes, kernels, and windows

LocalFilters.kernelFunction
kernel([Dims{N},] args...)

yields an N-dimensional abstract array built from args... and which can be used as a kernel in local filtering operations. The kernel is called a neighborhood, a sliding window, or a structuring element when its element are of type Bool and indicate whether a cell is to be considered. A box is a kernel whose elements are all true, it corresponds to an hyper-rectangular neighborhood whose sides are aligned with the Cartesian axes.

  • If args... is composed of N integers and/or ranges or if it is an N-tuple of integers and/or ranges, a box is returned whose axes are specified by args.... Each integer argument is converted in a centered unit range of this length (see LocalFilters.kernel_range).

  • If Dims{N} is provided and args... is a single integer or range, it is interpreted as being the same for all dimensions of an N-dimensional kernel. For example, kernel(Dims{3},5) yields a 3-dimensional box with index range -2:2 in every dimension.

  • If args... is two Cartesian indices or a 2-tuple of Cartesian indices, say start and stop, a box is returned whose first and last indices are start and stop.

  • If args... is a Cartesian range, say R::CartesianIndices{N}, a box array is returned whose axes are given by R.

  • If args... is an abstract array of any other type than CartesianIndices, it is returned unchanged.

Optional leading argument Dims{N} can be specified to assert the number of dimensions of the result or to provide the number of dimensions when it cannot be inferred from the arguments. For example, when args... is a single integer length or range which should be interpreted as being the same for all dimensions.

See also LocalFilters.strel, LocalFilters.ball, LocalFilters.kernel_range, and LocalFilters.reverse_kernel.

source
LocalFilters.KernelType
LocalFilters.Kernel{N}

is the union of types of arguments suitable to define a N-dimensional kernel or filter. Any argument B of this type can be converted into a N-dimensional kernel array by kernel(Dims{N}, B).

source
LocalFilters.WindowType
LocalFilters.Window{N}

is the union of types of arguments suitable to define a N-dimensional sliding window that is a N-dimensional array of Booleans indicating whether a position belongs to the window of not. Any argument B of this type can be converted into a N-dimensional array of Booleans by kernel(Dims{N}, B).

source
LocalFilters.strelFunction
strel(T, A)

yields a structuring element suitable for mathematical morphology operations. The result is an array whose elements have type T (which can be Bool or a floating-point type). Argument A can be a hyper-rectangular Cartesian sliding window or an array with boolean elements.

If T is a floating-point type, then the result is a so-called flat structuring element whose coefficients are zero(T) inside the shape defined by A and -T(Inf) elsewhere.

See also LocalFilters.kernel, LocalFilters.box, and LocalFilters.ball.

source
LocalFilters.ballFunction
LocalFilters.ball(Dims{N}, r)

yields a mask approximating a N-dimensional ball of radius r. The result is N-dimensional array of Boolean's with all dimensions odd and equal and whose values are true inside the ball (where distance to the center ≤ r) and false otherwise. The mask may be used to specify the neighborhood, the kernel, or the structuring element in local filtering operations.

The returned mask has centered axes, to get a mask with 1-based indices, call:

LocalFilters.ball(Dims{N}, r).parent

See also LocalFilters.kernel and LocalFilters.strel.

source
LocalFilters.BoxType
const LocalFilters.Box{N} = FastUniformArray{Bool,N,true}

is an alias to the type of uniformly true N-dimensional arrays. Instances of this kind are used to represent hyper-rectangular sliding windows in LocalFilters.

Method LocalFilters.box can be used to build an object of this type.

source
LocalFilters.centeredFunction

LocalFilters.centered(A) -> B

yields an abstract array B sharing the entries of array A but with offsets on indices so that the axes of B are centered (for even dimension lengths, the same conventions as in fftshift are used).

This public method is purposely not exported because it could introduce some confusions. For example OffsetArrays.centered is similar but has a slightly different semantic.

See also LocalFilters.kernel_range, LocalFilters.centered_offset.

source
LocalFilters.reverse_kernelFunction
B = reverse_kernel(args...)

yields a kernel B which is similar to A = kernel(args...) but with reversed ordering in the sense that B[i] == A[-i] holds for all indices i such that -i is a valid index in A. As a consequence, a correlation by B yields the same result as a convolution by A and conversely.

See also LocalFilters.kernel and LocalFilters.strel.

source

Utilities

Below are described non-exported types and methods that may be useful in the context of building local filters.

LocalFilters.IndicesType
I = LocalFilters.Indices(A::AbstractArray...)

builds a callable object, I, that can be used to produce ranges of indices for each of the arrays A.... These ranges will all be of the same type: linear index ranges, if all arrays A... are vectors, Cartesian index ranges otherwise.

I is similar to the eachindex method but is specialized for a style of indexing, it can be called as I(B...) to yield a suitable index range to access all the entries of array(s) B... which are any number of the A... specified when building I. If B... consists in several arrays, they must have the same axes otherwise I(B...) will throw a DimensionMismatch exception.

Call:

I = LocalFilters.Indices{S}()

with S = IndexLinear or S = IndexCartesian to specifically choose the indexing style.

source
LocalFilters.localindicesFunction
LocalFilters.localindices(A_inds, ord, B_inds, i) -> J

yields the subset J of all indices j such that:

  • A[j] and B[ord(i,j)] = B[j-i] are in-bounds if ord = FORWARD_FILTER;

  • A[j] and B[ord(i,j)] = B[i-j] are in-bounds if ord = REVERSE_FILTER;

with A and B arrays whose index ranges are given by A_inds and B_inds. To make the code agnostic to the ordering, use A[i] and B[ord(i,j)] to retrieve the values in A and B.

Index ranges A_inds and B_inds and index i must be of the same kind:

  • linear index ranges for A_inds and B_inds and linear index for i;

  • Cartesian index ranges for A_inds and B_inds and Cartesian index for i of same number of dimensions.

Constructor LocalFilters.Indices may by used to build a callable object that yields the index ranges of A and B in a consistent way:

indices = LocalFilters.Indices(A, B)
A_inds = indices(A)
B_inds = indices(B)
source
LocalFilters.check_axesFunction
LocalFilters.check_axes([I,] A...)

throws an exception if not all arrays A... have the same axes, or all have axes I if specified.

source
LocalFilters.is_morpho_math_boxFunction
LocalFilters.is_morpho_math_box(B)

yields whether structuring element B has the same effect as an hyper-rectangular box for mathematical morphology operations. This may be used to use fast separable versions of mathematical morphology operations like the van Herk-Gil-Werman algorithm.

source
LocalFilters.kernel_rangeFunction
LocalFilters.kernel_range([ord=FORWARD_FILTER,] rng::AbstractRange{<:Integer})
LocalFilters.kernel_range([ord=FORWARD_FILTER,] len::Integer)
LocalFilters.kernel_range([ord=FORWARD_FILTER,] start::Integer, stop::Integer)

yield an unit-step Int-valued index range based on range rng, dimension length len, or first and last indices start and stop. In the case of a given dimension length, a centered range of this length is returned (for even lengths, the same conventions as in fftshift are used).

If ordering ord is specified, the returned range is suitable for this ordering.

See also LocalFilters.kernel, LocalFilters.centered_offset, and LocalFilters.centered.

source
LocalFilters.unit_rangeFunction
LocalFilters.unit_range(r::Union{AbstractRange{<:Integer},CartesianIndices})

converts r into an Int-valued unit step index range. r may be a linear or a Cartesian index range. If r is a linear range, the absolute value of its step must be 1.

LocalFilters.unit_range(start::Integer, stop::Integer)

yields the Int-valued unit step range Int(start):Int(stop).

source