Package References

High-level API

ImageBase.restrictFunction
restrict(img[, dims]) -> imgr

Reduce the size of img by approximately two-fold along the dimensions listed in dims, or all spatial coordinates if dims is not specified.

Output

The type of output array imgr depends on the input type:

  • If img is an OffsetArray, then output array imgr will also be an OffsetArray.
  • If img is not an OffsetArray, then output array imgr will be an Array type even if it has offset indices.

The size of imgr is approximately 1/2 of the original size. More specifically:

  • if Nₖ = size(img, k) is odd, then size(imgr, k) == (Nₖ+1) ÷ 2.
  • if Nₖ = size(img, k) is even, then size(imgr, k) == (Nₖ÷2) + 1.

Examples

The optional argument dims can be a Tuple or Integer:

A = rand(5, 5) # size: (5, 5)

restrict(A) # size: (3, 3)

restrict(A, 1) # size: (3, 5)
restrict(A, 2) # size: (5, 3)

restrict(A, (1, )) # size: (3, 5)
restrict(A, (1, 2)) # size: (3, 3)

Unless the input array is 1-based, the origin will be halfed:

julia> using ImageBase, OffsetArrays

julia> Ao = OffsetArray(rand(5, 4), 5, 6);

julia> Ar = restrict(Ao);

julia> axes(Ao)
(OffsetArrays.IdOffsetRange(values=6:10, indices=6:10), OffsetArrays.IdOffsetRange(values=7:10, indices=7:10))

julia> axes(Ar)
(OffsetArrays.IdOffsetRange(values=3:5, indices=3:5), OffsetArrays.IdOffsetRange(values=4:6, indices=4:6))

Extended help

The term restrict is taken from the coarsening operation of algebraic multigrid methods; it is the adjoint of "prolongation" (which is essentially interpolation). restrict anti-aliases the image as it goes, so is better than a naive summation over 2x2 blocks. The implementation of restrict has been tuned for performance, and should be a fast method for constructing pyramids.

If l is the size of img along a particular dimension, restrict produces an array of size (l+1)÷2 for odd l, and l÷2 + 1 for even l. See the example below for an explanation.

See also ImageTransformations.imresize.

Example

a_course = [0, 1, 0.3]

If we were to interpolate this at the halfway points, we'd get

a_fine = [0, 0.5, 1, 0.65, 0.3]

Note that a_fine is obtained from a_course via the prolongation operator P as P*a_course, where

P = [1   0   0;      # this line "copies over" the first point
     0.5 0.5 0;      # this line takes the mean of the first and second point
     0   1   0;      # copy the second point
     0   0.5 0.5;    # take the mean of the second and third
     0   0   1]      # copy the third

restrict is the adjoint of prolongation. Consequently,

julia> restrict(a_fine)
3-element Array{Float64,1}:
 0.125
 0.7875
 0.3125

julia> (P'*a_fine)/2
3-element Array{Float64,1}:
 0.125
 0.7875
 0.3125

where the division by 2 approximately preserves the mean intensity of the input.

As we see here, for odd-length a_fine, restriction is the adjoint of interpolation at half-grid points. When length(a_fine) is even, restriction is the adjoint of interpolation at 1/4 and 3/4-grid points. This turns out to be the origin of the l->l÷2 + 1 behavior.

One consequence of this definition is that the edges move towards zero:

julia> restrict(ones(11))
6-element Array{Float64,1}:
 0.75
 1.0
 1.0
 1.0
 1.0
 0.75

In some applications (e.g., image registration), you may find it useful to trim the edges.

ImageTransformations.imrotateFunction
imrotate(img, θ, [indices]; kwargs...) -> imgr

Rotate image img by θ∈[0,2π) in a clockwise direction around its center point.

Arguments

  • img::AbstractArray: the original image that you need to rotate.
  • θ::Real: the rotation angle in clockwise direction. To rotate the image in conter-clockwise direction, use a negative value instead. To rotate the image by d degree, use the formular θ=d*π/180.
  • indices (Optional): specifies the output image axes. By default, rotated image imgr will not be cropped, and thus axes(imgr) == axes(img) does not hold in general.

Parameters

Info

To construct method and fillvalue values, you may need to load Interpolations package first.

  • method::Union{Degree, InterpolationType}: the interpolation method you want to use. By default it is Linear().
  • fillvalue: the value that used to fill the new region. The default value is NaN if possible, otherwise is 0.

This function is a simple high-level interface to warp, for more explaination and details, please refer to warp.

Examples

using TestImages, ImageTransformations
img = testimage("cameraman")

# Rotate the image by π/4 in the clockwise direction
imgr = imrotate(img, π/4) # output axes (-105:618, -105:618)

# Rotate the image by π/4 in the counter-clockwise direction
imgr = imrotate(img, -π/4) # output axes (-105:618, -105:618)

# Preserve the original axes
# Note that this is more efficient than `@view imrotate(img, π/4)[axes(img)...]`
imgr = imrotate(img, π/4, axes(img)) # output axes (1:512, 1:512)

By default, imrotate uses bilinear interpolation with constant fill value (NaN or 0). You can, for example, use the nearest interpolation and fill the new region with white pixels:

using Interpolations, ImageCore
imrotate(img, π/4, method=Constant(), fillvalue=oneunit(eltype(img)))

And with some inspiration, maybe fill with periodic values and tile the output together to get a mosaic:

using Interpolations, ImageCore
imgr = imrotate(img, π/4, fillvalue = Periodic())
mosaicview([imgr for _ in 1:9]; nrow=3)
source
ImageTransformations.imresizeFunction
imresize(img, sz; [method]) -> imgr
imresize(img, inds; [method]) -> imgr
imresize(img; ratio, [method]) -> imgr

upsample/downsample the image img to a given size sz or axes inds using interpolations. If ratio is provided, the output size is then ceil(Int, size(img).*ratio).

Tip

This interpolates the values at sub-pixel locations. If you are shrinking the image, you risk aliasing unless you low-pass filter img first.

Arguments

  • img: the input image array
  • sz: the size of output array
  • inds: the axes of output array If inds is passed, the output array imgr will be OffsetArray.

Parameters

Info

To construct method, you may need to load Interpolations package first.

  • ratio: the upsample/downsample ratio used. The output size is ceil(Int, size(img).*ratio). If ratio is larger than 1, it is an upsample operation. Otherwise it is a downsample operation. ratio can also be a tuple, in which case ratio[i] specifies the resize ratio at dimension i.
  • method::InterpolationType: specify the interpolation method used for reconstruction. conveniently, methold can also be a Degree type, in which case a BSpline object will be created. For example, method = Linear() is equivalent to method = BSpline(Linear()).

Examples

using ImageTransformations, TestImages, Interpolations

img = testimage("lighthouse") # 512*768

# pass integers as size
imresize(img, 256, 384) # 256*384
imresize(img, (256, 384)) # 256*384
imresize(img, 256) # 256*768

# pass indices as axes
imresize(img, 1:256, 1:384) # 256*384
imresize(img, (1:256, 1:384)) # 256*384
imresize(img, (1:256, )) # 256*768

# pass resize ratio
imresize(img, ratio = 0.5) #256*384
imresize(img, ratio = (2, 1)) # 1024*768

# use different interpolation method
imresize(img, (256, 384), method=Linear()) # 256*384 bilinear interpolation
imresize(img, (256, 384), method=Lanczos4OpenCV()) # 256*384 OpenCV-compatible Lanczos 4 interpolation

For downsample with ratio=0.5, restrict is a much faster two-fold implementation that you can use.

source

Low-level warping API

ImageTransformations.warpFunction
warp(img, tform, [indices]; kwargs...) -> imgw

Transform the coordinates of img, returning a new imgw satisfying imgw[I] = img[tform(I)].

Output

The output array imgw is an OffsetArray. Unless manually specified, axes(imgw) == axes(img) does not hold in general. If you just want a plain array, you can "strip" the custom indices with parent(imgw) or OffsetArrays.no_offset_view(imgw).

Arguments

  • img: the original image that you need coordinate transformation.
  • tform: the coordinate transformation function or function-like object, it must accept a SVector as input. A useful package to create a wide variety of such transfomrations is CoordinateTransformations.jl.
  • indices (Optional): specifies the output image axes. By default, the indices are computed in such a way that imgw contains all the original pixels in img using autorange. To do this inv(tform) has to be computed. If the given transfomration tform does not support inv then the parameter indices has to be specified manually.

Parameters

Info

To construct method and fillvalue values, you may need to load Interpolations package first.

  • method::Union{Degree, InterpolationType}: the interpolation method you want to use. By default it is BSpline(Linear()). To construct the method instance, one may need to load Interpolations.
  • fillvalue: the value that used to fill the new region. The default value is NaN if possible, otherwise is 0. One can also pass the extrapolation boundary condition: Flat(), Reflect() and Periodic().

See also

There're some high-level interfaces of warp:

There are also lazy versions of warp:

Extended help

Parameters in detail

This approach is known as backward mode warping. It is called "backward" because the internal coordinate transformation is actually an inverse map from axes(imgr) to axes(img).

You can manually specify interpolation behavior by constructing AbstractExtrapolation object and passing it to warp as img. However, this is usually cumbersome. For this reason, there are two keywords method and fillvalue to conveniently construct an AbstractExtrapolation object during warp.

Warning

If img is an AbstractExtrapolation, then additional method and fillvalue keywords will be discarded.

method::Union{Degree, InterpolationType}

The interpolation method you want to use to reconstruct values in the wrapped image.

Among those possible InterpolationType choice, there are some commonly used methods that you may have used in other languages:

  • nearest neighbor: BSpline(Constant())
  • triangle/bilinear: BSpline(Linear())
  • bicubic: BSpline(Cubic(Line(OnGrid())))
  • lanczos2: Lanczos(2)
  • lanczos3: Lanczos(3)
  • lanczos4: Lanczos(4) or Lanczos4OpenCV()

When passing a Degree, it is expected to be a BSpline. For example, Linear() is equivalent to BSpline(Linear()).

fillvalue

In case tform(I) maps to indices outside the original img, those locations are set to a value fillvalue. The default fillvalue is NaN if the element type of img supports it, and 0 otherwise.

The parameter fillvalue can be either a Number or Colorant. In this case, it will be converted to eltype(imgr) first. For example, fillvalue = 1 will be converted to Gray(1) which will fill the outside indices with white pixels.

Also, fillvalue can be extrapolation schemes: Flat(), Periodic() and Reflect(). The best way to understand these schemes is perhaps try it with small example:

using ImageTransformations, TestImages, Interpolations
using OffsetArrays: IdOffsetRange

img = testimage("lighthouse")

imgr = imrotate(img, π/4; fillvalue=Flat()) # zero extrapolation slope
imgr = imrotate(img, π/4; fillvalue=Periodic()) # periodic boundary
imgr = imrotate(img, π/4; fillvalue=Reflect()) # mirror boundary

axes(imgr)

# output

(IdOffsetRange(values=-196:709, indices=-196:709), IdOffsetRange(values=-68:837, indices=-68:837))

The meaning of the coordinates

imgw keeps track of the indices that would result from applying inv(tform) to the indices of img. This can be very handy for keeping track of how pixels in imgw line up with pixels in img.

using ImageTransformations, TestImages, Interpolations

img = testimage("lighthouse")
imgr = imrotate(img, π/4)
imgr_cropped = imrotate(img, π/4, axes(img))

# No need to manually calculate the offsets
imgr[axes(img)...] == imgr_cropped

# output
true
Tip

For performance consideration, it's recommended to pass the inds positional argument to warp instead of cropping the output with imgw[inds...].

Examples: a 2d rotation

Note

This example only shows how to construct tform and calls warp. For common usage, it is recommended to use imrotate function directly.

Rotate around the center of img:

using ImageTransformations, CoordinateTransformations, Rotations, TestImages, OffsetArrays
using OffsetArrays: IdOffsetRange
img = testimage("lighthouse") # axes (1:512, 1:768)

tfm = recenter(RotMatrix(-pi/4), center(img))
imgw = warp(img, tfm)

axes(imgw)

# output

(IdOffsetRange(values=-196:709, indices=-196:709), IdOffsetRange(values=-68:837, indices=-68:837))
source
ImageTransformations.WarpedViewType
WarpedView(img, tform, [indices]; kwargs...) -> wv

Create a view of img that lazily transforms any given index I passed to wv[I] so that wv[I] == img[tform(I)].

This is the lazy view version of warp, please see warp for more information.

source
ImageTransformations.InvWarpedViewType
InvWarpedView(img, tinv, [indices]; kwargs...) -> wv
InvWarpedView(inner_view, tinv) -> wv

Create a view of img that lazily transforms any given index I passed to wv[I] so that wv[I] == img[inv(tinv)(I)].

Except for the lazy evaluation, the following two lines are expected to be equivalent:

warp(img, inv(tform), [indices]; kwargs...)
invwarpedview(img, tform, [indices]; kwargs...)

The conceptual difference to WarpedView is that InvWarpedView is intended to be used when reasoning about the image is more convenient that reasoning about the indices. Furthermore, InvWarpedView allows simple nesting of transformations, in which case the transformations will be composed into a single one.

For detailed explaination of warp, associated arguments and parameters, please refer to warp.

source

Utilities

ImageTransformations.autorangeFunction
autorange(A::AbstractArray, tform::Transformation)

For given transformation tform, return the "smallest" range indices that preserves all information from A after applying tform.

Examples

For transformation that preserves the array size, autorange is equivalent to axes(A).

A = rand(5, 5)
tform = IdentityTransformation()
autorange(A, tform) == axes(A)

# output
true

The diffrence shows up when tform enlarges the input array A. In the following example, we need at least (0:6, 0:6) as the range indices to get all data of A:

A = rand(5, 5)
tform = recenter(RotMatrix(pi/8), center(A))
autorange(A, tform)

# output
(0:6, 0:6)
Note

This function is not exported; it is mainly for internal usage to infer the default indices.

source