Contents
Vigranumpy exports the functionality of the C++ image processing library VIGRA to Python. It can be invoked by importing the vigra module:
import vigra
Vigranumpy is based on the popular numpy module and uses its ndarray data structure to represent image and volume data. It introduces the C++ wrapper class NumpyArray to allow transparent execution of VIGRA C++ functions on numpy arrays. Thus, vigranumpy is fully interoperable with existing numpy functionality, including various tools for image display such as matplotlib. Since vigranumpy uses boost_python, it is able to support function overloading (which plain Python does not provide), so that calling syntax is largely uniform, regardless of the type and dimension of the input arguments.
Basic calling syntax is similar to C++, with one important difference: Arguments for output images are optional. If no output image is provided, vigranumpy will allocate it as appropriate. In either case, the output image will be returned by the function, for example:
# allocate new result image
>>> smoothImage = vigra.gaussianSmoothing(inputImage, scale)
# reuse and overwrite existing result image
>>> smoothImage = vigra.gaussianSmoothing(inputImage, scale, out=smoothImage)
Unless otherwise noted, all functions expect and create arrays with dtype=numpy.float32.
Another important concern is the interpretation and ordering of the array axes. Numpy does not provide any means to attach semantics to axes, but relies purely on the convention that the most important axis is last, as in array[y, x] or array[z, y, x] (“C-order”). However, there is no way to enforce this convention in a program, since arrays can be transposed outside of the user’s control (e.g. when saving data to a file). Moreover, many imaging libraries (e.g. Image Magick, OpenCV, Qt and the C++ version of VIGRA) use the opposite convention where the x-axis comes first, as in array[x, y] or array[x, y, z]. This makes it very difficult to convert solutions developed in Python into a fast C++ version, because one has to reverse all indices without making mistakes. Matters become even more complicated when multi-channel (e.g. color) images are considered – should the color index now be first or last?
To solve these ambiguities in a clean way, vigranumpy introduces the concept of axistags which is realized in class vigra.AxisTags. Every VigraArray (which is a subclass of numpy.ndarray) gets a new property array.axistags that describes axis semantics, and all vigranumpy functions account for and preserve axistag information. Unfortunately, this functionality cannot easily be retrofitted to numpy.ndarray itself. Therefore, we employ the following conversion rules between Python and C++ arrays:
When the Python array has no array.axistags property, it is mapped to the C++ NumpyArray without any change in axis ordering. Since most VIGRA functions can work on arbitrarily transposed arrays, you will get correct results, but execution may be slower because the processor cache is poorly utilized in certain axis orders.
Moreover, this may lead to overload resolution ambiguities. For example, when the array has shape (3, 60, 40), vigranumpy has no way to decide if this is a 2-dimensional RGB image or a 3-dimensional array that happens to have only 3 slices. Thus, vigranumpy may not always execute the function you actually intended to call.
When the Python array has the array.axistags property, it is transposed into a canonical axis ordering before vigranumpy executes a function, and the results are transposed back into the original ordering. Likewise, functions that change axis ordering (such as array.swapaxes(0,1)) or reduce the number of axes (such as array.max(axis=1)) as well as array arithmetic operations preserve axistags (see section Mathematical Functions and Type Coercion). Thus, you can work in any desired axis order without loosing control. Overload ambiguities can no longer occur because a function cannot be called when the axistags are unsuitable.
Detailed information about the use of axistags is given in section Axistags and the VigraArray Data Structure below. Section Writing Your Own C++ Modules describes how you can take advantage of the axistags mechanism in your own C++ code.
While vigranumpy can directly work on numpy.ndarrays, this would not give us the advantages of axistags as described above. Therefore, vigranumpy introduces its own array class VigraArray which is a subclass of numpy.ndarray, but re-implements many of its methods so that axistags are respected. Arrays with a conforming axistags property are most easily constructed by one of the predefined array factories. We illustrate the ideas by some examples:
>>> width, height, depth = 300, 200, 3
# create a 2-dimensional RGB image
>>> rgb = vigra.RGBImage((width, height))
>>> rgb.shape
(300, 200, 3)
>>> rgb.axistags # short output: only axis keys
x y c
>>> print rgb.axistags # long output
AxisInfo: 'x' (type: Space)
AxisInfo: 'y' (type: Space)
AxisInfo: 'c' (type: Channels) RGB
# create a 3-dimensional scalar volume
>>> volume = vigra.ScalarVolume((width, height, depth))
>>> volume.shape
(300, 200, 3) # same shape as before
>>> volume.axistags
x y z # but different semantic interpretation
>>> print volume.axistags
AxisInfo: 'x' (type: Space)
AxisInfo: 'y' (type: Space)
AxisInfo: 'z' (type: Space)
It is also possible to attach additional information to the axistags, in particular the resolution of the axis, and a text comment. The resolution will be correctly adjusted when the image is resized:
>>> rgb.axistags['x'].resolution = 1.2 # in some unit of length
>>> rgb.axistags['y'].resolution = 1.4 # in some unit of length
>>> rgb.axistags['c'].description = 'fluorescence microscopy, DAPI and GFP staining'
>>> print rgb.axistags
AxisInfo: 'x' (type: Space, resolution=1.2)
AxisInfo: 'y' (type: Space, resolution=1.4)
AxisInfo: 'c' (type: Channels) fluorescence microscopy, DAPI and GFP staining
# interpolate the image to twice its original size
>>> rgb2 = vigra.sampling.resize(rgb, shape=(2*width-1, 2*height-1))
>>> print rgb2.axistags
AxisInfo: 'x' (type: Space, resolution=0.6)
AxisInfo: 'y' (type: Space, resolution=0.7)
AxisInfo: 'c' (type: Channels) fluorescence microscopy, DAPI and GFP staining
When the array is transposed, the axistags are transposed accordingly. When axes are dropped from the array, the corresponding entries are dropped from the axistags property:
# transpose the volume into reverse axis order
>>> transposed_volume = volume.transpose()
>>> transposed_volume.shape
(3, 200, 300)
>>> transposed_volume.axistags
z y x
# get a view to the first slice (z == 0)
>>> first_slice = volume[..., 0]
>>> first_slice.shape
(300, 200)
>>> first_slice.axistags
x y
# get the maximum of each slice
>>> volume.max(axis=0).max(axis=0)
VigraArray(shape=(3,), axistags=z, dtype=float32, data=
[ 0. 0. 0.])
# likewise, but specify axes by their keys
>>> volume.max(axis='x').max(axis='y')
VigraArray(shape=(3,), axistags=z, dtype=float32, data=
[ 0. 0. 0.])
The initial ordering of the axes is controlled by the argument order that can optionally be passed to the VigraArray constuctor or the factory functions. If order is not explicitly provided, it is determined by the static property VigraArray.defaultOrder (which yields ‘V’ order). The following orders are currently supported:
- ‘C’ order:
- Both strides and axes are arranged in descending order, as in a plain numpy.ndarray. For example, axistags will be ‘y x c’ or ‘z y x c’. array.flags[‘C_CONTIGUOUS’] will be true.
- ‘F’ order:
- Both strides and axes are arranged in ascending order, i.e. opposite to ‘C’ order. For example, axistags will be ‘c x y’ or ‘c x y z’. array.flags[‘F_CONTIGUOUS’] will be true.
- ‘V’ order:
- VIGRA-order is an interleaved memory layout that simulates vector-valued pixels or voxels: Channels will be the last axis and have the smallest stride, whereas all other axes are arranged in ascending order. For example, axistags will be ‘x y c’ or ‘x y z c’.
- ‘A’ order:
- Defaults to ‘V’ when a new array is created, and means ‘preserve order’ when an existing array is copied.
The meaning of ‘ascending’ or ‘descending’ order is determined by two rules: the primary order is according to axis type (see vigra.AxisType), where Channels < Space < Angle < Time < Frequency < Unknown. The secondary order (between axes of the same type) is lexicographic, such that ‘x’ < ‘y’ < ‘z’. Usage examples:
>>> rgbv = vigra.RGBImage((width, height), order='V')
>>> rgbv.shape
(300, 200, 3)
>>> rgbv.axistags
x y c
>>> rgbc = vigra.RGBImage((width, height), order='C')
>>> rgbc.shape
(200, 300, 3)
>>> rgbc.axistags
y x c
>>> rgbf = vigra.RGBImage((width, height), order='F')
>>> rgbf.shape
(3, 300, 200)
>>> rgbf.axistags
c x y
Functions that reduce the array to a one-dimensional shape (flatten(), flat, ravel(), take()) always transpose the array into ‘C’ order before flattening.
Axistags are stored in a list-like class vigra.AxisTags, whose individual entries are of type vigra.AxisInfo. The simplest way to attach axistags to a plain numpy.ndarray (by creating a view of type VigraArray) is via the convenience function vigra.taggedView().
Enum to encode the type of an axis described by an AxisInfo object. Possible values:
- AxisType.Channels:
- a channel axis
- AxisType.Space:
- a spatial axis
- AxisType.Angle:
- an axis encoding angles (e.g. polar coordinates)
- AxisType.Time:
- a temporal axis
- AxisType.Frequency:
- an axis in the Fourier domain
- AxisType.UnknownAxisType:
- type not specified
- AxisType.NonChannel:
- any type except Channels
- AxisType.AllAxes:
- any type
Types can be combined by the bitwise ‘or’ operator. For example, Space | Frequency denotes a spatial axis in the Fourier domain.
An object describing a single axis.
Constructor:
Parameters: |
|
---|
In addition, AxisInfo defines the following factories for the most common cases:
- AxisInfo.c or AxisInfo.c(description='a description'):
- Factory for an axisinfo object describing the ‘c’ (channel) axis.
- AxisInfo.x or AxisInfo.x(resolution=0.0, description=''):
- Factory for an axisinfo object describing the ‘x’ (spatial) axis.
- AxisInfo.y or AxisInfo.y(resolution=0.0, description=''):
- Factory for an axisinfo object describing the ‘y’ (spatial) axis.
- AxisInfo.z or AxisInfo.z(resolution=0.0, description=''):
- Factory for an axisinfo object describing the ‘z’ (spatial) axis.
- AxisInfo.t or AxisInfo.t(resolution=0.0, description=''):
- Factory for an axisinfo object describing the ‘t’ (time) axis.
- AxisInfo.fx or AxisInfo.fx(resolution=0.0, description=''):
- Factory for an axisinfo object describing the ‘x’ axis in the Fourier domain.
- AxisInfo.fy or AxisInfo.fy(resolution=0.0, description=''):
- Factory for an axisinfo object describing the ‘y’ axis in the Fourier domain.
- AxisInfo.fz or AxisInfo.fz(resolution=0.0, description=''):
- Factory for an axisinfo object describing the ‘z’ axis in the Fourier domain.
- AxisInfo.ft or AxisInfo.ft(resolution=0.0, description=''):
- Factory for an axisinfo object describing the ‘t’ axis in the Fourier domain.
(read/write property, type ‘float’) the resolution of the axis. The resolution will be automatically adjusted whenever the image size changes, e.g. due to a call to resize() or slicing with non-unit step size:
>>> a = vigra.RGBImage((200,100))
>>> a.axistags['x'].resolution = 1.0
>>> a.axistags['y'].resolution = 1.2
>>> print a.axistags
AxisInfo: 'x' (type: Space, resolution=1)
AxisInfo: 'y' (type: Space, resolution=1.2)
AxisInfo: 'c' (type: Channels) RGB
>>> b = a[::2, ::4, :]
>>> print b.axistags
AxisInfo: 'x' (type: Space, resolution=2)
AxisInfo: 'y' (type: Space, resolution=4.8)
AxisInfo: 'c' (type: Channels) RGB
Object to describe axis properties and axis ordering in a VigraArray.
Constructor:
Most AxisTags methods should not be called directly, but via the corresponding array methods, because this ensures that arrays and their axistags are always kept in sync (rule of thumb: if an axistags function is not documented, you call it on your own risk).
The entries of an axistags object (i.e. the individual axisinfo objects) can be accessed via the index operator, where the argument can either be the axis index or the axis key:
>>> print array.axistags[0]
AxisInfo: 'x' (type: Space, resolution=1.2)
>>> print array.axistags['x']
AxisInfo: 'x' (type: Space, resolution=1.2)
>>> array.axistags['x'].resolution = 2.0
>>> print array.axistags['x']
AxisInfo: 'x' (type: Space, resolution=2)
Get the axis index of a given axis key:
>>> axistags.index('x')
0 # the 'x' axis is first index
How many axes of the given type(s) are in this axistags object?:
axistags.axisTypeCount(types) -> int
The ‘types’ of the query must be single AxisType instances or a combination of them. Examples:
>>> a = vigra.defaultAxistags('txyc')
>>> a.axisTypeCount(vigra.AxisType.Space)
2
>>> a.axisTypeCount(vigra.AxisType.Time)
1
>>> a.axisTypeCount(vigra.AxisType(vigra.AxisType.Space | vigra.AxisType.Time))
3
>>> a.axisTypeCount(vigra.AxisType.NonChannel)
3
Set a description for the channel axis, if one exists:
axistags.setChannelDescription('colors are in Lab color space')
It is similar to:
axistags['c'].description = 'colors are in Lab color space'
except when the axistags contain no channel axis, in which case setChannelDescription() is simply ignored, whereas axistags[‘c’] would cause an exception.
Bases: numpy.ndarray
This class extends numpy.ndarray with the concept of axistags which encode the semantics of the array’s axes. VigraArray overrides all numpy.ndarray methods in order to handle axistags in a sensible way. In particular, operations acting on two arrays simultaneously (e.g. addition) will first transpose the arguments such that their axis ordering matches.
Constructor:
Parameters: |
|
---|
obj may be one of the following
order can be ‘C’ (C order), ‘F’ (Fortran order), ‘V’ (VIGRA order), ‘A’ (any), or None. This parameter controls the order of strides and axistags (unless axistags are explicit passed into the constructor). See the order definitions for details. If ‘order=None’, the order is determined by VigraArray.defaultOrder.
Get default axistags for the given specification ‘tagSpec’. TagSpec can be the number of dimensions of the array (array.ndim, must be <= 5) or a string containing a sequence of axis keys (only the default keys ‘x’, ‘y’, ‘z’, ‘t’, and ‘c’ are currently supported). The ‘order’ parameter determines the axis ordering, see the order definitions for details. If ‘noChannels’ is True, there will be no channel axis. Examples:
>>> vigra.VigraArray.defaultAxistags(3)
x y c
>>> vigra.VigraArray.defaultAxistags(4)
x y z c
>>> vigra.VigraArray.defaultAxistags(5)
x y z t c
>>> vigra.VigraArray.defaultAxistags(3, order='C')
y x c
>>> vigra.VigraArray.defaultAxistags(2, noChannels=True)
x y z
>>> vigra.VigraArray.defaultAxistags(3, noChannels=True)
x y z
>>> vigra.VigraArray.defaultAxistags(4, noChannels=True)
x y z t
>>> vigra.VigraArray.defaultAxistags('xty')
x t y
>>> vigra.VigraArray.defaultAxistags('xty', order='V')
x y t
Create a view containing the desired axis keys in the given order. When the array contains an axis not listed, the axis will be dropped if it is a singfleton (otherwise, an exception is raised). If a requested key is not present in this array, a singleton axis will be inserted at that position, if the missing key is among the known standard keys (otherwise, an exception is raised). The function fails if this array contains axes of unknown type (key ‘?’). If ‘self’ is already suitable, it is simply retured without generating a new view.
Usage:
>>> a = vigra.ScalarVolume((200, 100))
>>> a.axistags
x y
>>> a.shape
(200, 100)
>>> b = a.withAxes('y', 'x', 'c')
>>> b.axistags
y x c
>>> b.shape
(100, 200, 1)
array.__getitem__(index) implements the indexing operator array[index]. In addition to the usual numpy.ndarray indexing functionality, this function also updates the axistags of the result array. There are three cases:
- getitem creates a scalar value => no axistags are required
- getitem creates an array view => axistags are transferred from the corresponding axes of the base array
- getitem creates a copy of an array (fancy indexing) => all axistags are ‘?’
If the index contains ‘numpy.newaxis’, a new singleton axis is inserted at the appropriate position, whose axisinfo is set to ‘?’ (unknown). If the index contains ‘vigra.newaxis(axisinfo)’, the singleton axis will get the given axisinfo.
Bind the axis identified by ‘which’ to the given ‘index’. This is similar to:
array[:, index, ...]
but you don not need to know the position of the axis when you use the axis key (according to axistags). For example, to get the green channel of an RGBImage, you write:
>>> rgb = vigra.RGBImage((200, 100))
>>> green = rgb.bindAxis('c', 1)
This gives the correct result irrespective of the axis ordering.
Create an iterator over the channels of the array. In each iteration, you get the array corresponding to a single channel. If the axistags contain no channel axis, there is only one iteration which yields the entire array. Example:
>>> rgb = vigra.RGBImage((200, 100))
>>> rgb.axistags
x y c
>>> red, green, blue = rgb.channelIter()
>>> red.axistags
x y
>>> red.shape
(200, 100)
Create an iterator over a single spatial axis of the array. In each iteration, you get the array corresponding to one coordinate along the axis given by ‘key’. For example, to iterate along the z-axis to get all x-y-slices in turn, you write:
>>> volume = vigra.Volume((width, height, depth))
>>> for slice in volume.sliceIter('z'):
... processSlice(slice)
Create an iterator over all the spatial coordinates in the array. In each iteration, you get the value corresponding to a single coordinate location. If the axistags contain no spatial axes, there is only one iteration which yields the entire array. Example:
>>> s = vigra.ScalarImage((2,2))
>>> s.ravel()[...] = range(4)
>>> for p in s.spaceIter():
.... print p
0.0
1.0
2.0
3.0
Create an iterator over the time points of the array. In each iteration, you get the array corresponding to a single time point. If the axistags contain no time axis, there is only one iteration which yields the entire array. Example:
>>> from vigra import *
>>> axistags = AxisTags(AxisInfo.t, AxisInfo.x, AxisInfo.y)
>>> timesteps, width, height = 2, 200, 100
>>> image_sequence = Image((timesteps, width, height), axistags=axistags)
>>> step1, step2 = image_sequence.timeIter()
Copy the values of an array to another one. This is similar to:
self[...] = other
but will first transpose both arrays so that axistags are aligned. If there is no valid alignment, RuntimeError will be raised.
Get a transposed view onto this array according to the given ‘order’. Possible orders are:
Convert this image to a Qt QImage (mainly for display purposes). The present image must have 1, 2, 3, or 4 channels, and the resulting QImage will have QImage.Format_Indexed8 iff there was only one channel and QImage.Format_[A]RGB32 otherwise (with the last of 2/4 channels being used as alpha channel).
The parameter normalize can be used to normalize an image’s value range to 0..255:
Display this image in a vigra.pyqt.ImageWindow.
The channels are intepreted as follows: 1 channel = gray image, 2 channels = gray + alpha, 3 channels = RGB, 4 channels = RGB + alpha.
The parameter normalize can be used to normalize an image’s value range to 0..255:
Create a new singleton axis via the indexing operator. This works similar to numpy.newaxis, but allows to provide an AxisInfo object for the new axis. For example:
>>> s = vigra.ScalarImage((width, height))
>>> s.axistags # no channel axis
x y
>>> t = s[..., numpy.newaxis]
>>> t.axistags # with unknown axis type
x y ?
>>> t = s[..., vigra.newaxis(vigra.AxisInfo.c)]
>>> t.axistags # with channel axis
x y c
Create a view to the given array with type VigraArray and the given axistags. This is a shorthand for:
>>> view = array.view(vigra.VigraArray)
>>> view.axistags = copy.copy(axistags)
The module vigra.impex defines read and write functions for image and volume data. Note that the contents of this module are automatically imported into the vigra module, so you may call ‘vigra.readImage(...)’ instead of ‘vigra.impex.readImage(...)’ etc.
Check whether the given file name contains image data:
isImage(filename) -> bool
This function tests whether a file has a supported image format. It checks the first few bytes of the file and compares them with the “magic strings” of each recognized image format. If the image format is supported it returns True otherwise False.
Ask for the image file extensions that vigra.impex understands:
listExtensions() -> string
This function returns a string containing the supported image file extensions for reading and writing with the functions readImage() and writeImage().
Ask for the image file formats that vigra.impex understands:
listFormats() -> string
This function returns a string containing the supported image file formats for reading and writing with the functions readImage() and writeImage().
Check how many images the given file contains:
numberImages(filename) -> int
This function tests how many images an image file contains(Values > 1 are only expected for the TIFF format to support multi-image TIFF).
Read an array from an HDF5 file.
‘filenameOrGroup’ can contain a filename or a group object referring to an already open HDF5 file. ‘pathInFile’ is the name of the dataset to be read, including intermediate groups. If the first argument is a group object, the path is relative to this group, otherwise it is relative to the file’s root group.
If the dataset has an attribute ‘axistags’, the returned array will have type VigraArray and will be transposed into the given ‘order’ (‘vigra.VigraArray.defaultOrder’ will be used if no order is given). Otherwise, the returned array is a plain ‘numpy.ndarray’. In this case, order=’F’ will return the array transposed into Fortran order.
Requirements: the ‘h5py’ module must be installed.
Read an image from a file:
readImage(filename, dtype = 'FLOAT', index = 0, order='') -> Image
When import_type is ‘UINT8’, ‘INT16’, ‘UINT16’, ‘INT32’, ‘UINT32’, ‘FLOAT’, ‘DOUBLE’, or one of the corresponding numpy dtypes (numpy.uint8 etc.), the returned image will have the requested pixel type. If dtype is ‘NATIVE’ or ‘’ (empty string), the image is imported with the original type of the data in the file. Caution: If the requested dtype is smaller than the original type in the file, values will be clipped at the bounds of the representable range, which may not be the desired behavior.
Individual images of sequential formats such as multi-image TIFF can be accessed via index. The number of images in a file can be checked with the function :func:`numberImages`(filename).
The order parameter determines the axis ordering of the resulting array (allowed values: ‘C’, ‘F’, ‘V’). When order == ‘’ (the default), vigra.VigraArray.defaultOrder is used.
Supported file formats are listed by the function listFormats(). When filename does not refer to a recognized image file format, an exception is raised. The file can be checked beforehand with the function isImage().
Read a 3D volume from a directory:
readVolume(filename, dtype='FLOAT', order='') -> Volume
If the volume is stored in a by-slice manner (e.g. one image per slice), the ‘filename’ can refer to an arbitrary image from the set. readVolume() then assumes that the slices are enumerated like:
name_base+[0-9]+name_ext
where name_base, the index, and name_ext are determined automatically. All slice files with the same name base and extension are considered part of the same volume. Slice numbers must be non-negative, but can otherwise start anywhere and need not be successive. Slices will be read in ascending numerical (not lexicographic) order. All slices must have the same size.
Otherwise, readVolume() will try to read ‘filename’ as an info text file with the following key-value pairs:
name = [short descriptive name of the volume] (optional)
filename = [absolute or relative path to raw voxel data file] (required)
gradfile = [abs. or rel. path to gradient data file] (currently ignored)
description = [arbitrary description of the data set] (optional)
width = [positive integer] (required)
height = [positive integer] (required)
depth = [positive integer] (required)
datatype = [UNSIGNED_CHAR | UNSIGNED_BYTE] (default: UNSIGNED_CHAR)
Lines starting with # are ignored. When import_type is ‘UINT8’, ‘INT16’, ‘UINT16’, ‘INT32’, ‘UINT32’, ‘FLOAT’, ‘DOUBLE’, or one of the corresponding numpy dtypes (numpy.uint8 etc.), the returned volume will have the requested pixel type.
The order parameter determines the axis ordering of the resulting array (allowed values: ‘C’, ‘F’, ‘V’). When order == ‘’ (the default), vigra.VigraArray.defaultOrder is used.
For details see the help for readImage().
Write an array to an HDF5 file.
‘filenameOrGroup’ can contain a filename or a group object referring to an already open HDF5 file. ‘pathInFile’ is the name of the dataset to be written, including intermediate groups. If the first argument is a group object, the path is relative to this group, otherwise it is relative to the file’s root group. If the dataset already exists, it will be replaced without warning.
If ‘data’ has an attribute ‘axistags’, the array is transposed to numpy order before writing. Moreover, the axistags will be stored along with the data in an attribute ‘axistags’.
Requirements: the ‘h5py’ module must be installed.
Save an image to a file:
writeImage(image, filename, dtype = '', compression = '', mode = 'w')
Parameters:
- image:
- the image to be saved
- filename:
- the file name to save to. The file type will be deduced from the file name extension (see vigra.impexListExtensions() for a list of supported extensions).
- dtype:
the pixel type written to the file. Possible values:
- ‘’ or ‘NATIVE’:
- save with original pixel type, or convert automatically when this type is unsupported by the target file format
- ‘UINT8’, ‘INT16’, ‘UINT16’, ‘INT32’, ‘UINT32’, ‘FLOAT’, ‘DOUBLE’:
- save as specified, or raise exception when this type is not supported by the target file format (see list below)
- ‘NBYTE’:
- normalize to range 0...255 and then save as ‘UINT8’
- numpy.uint8, numpy.int16 etc.:
- behaves like the corresponding string argument
- compression:
how to compress the data (ignored when compression type is unsupported by the file format). Possible values:
- ‘’ or not given:
- save with the native compression of the target file format
- ‘RLE’, ‘RunLength’:
- use run length encoding (native in BMP, supported by TIFF)
- ‘DEFLATE’:
- use deflate encoding (only supported by TIFF)
- ‘LZW’:
- use LZW algorithm (only supported by TIFF with LZW enabled)
- ‘ASCII’:
- write as ASCII rather than binary file (only supported by PNM)
- ‘1’ ... ‘100’:
- use this JPEG compression level (only supported by JPEG and TIFF)
- mode:
support for sequential file formats such as multi-image TIFF. Possible values:
‘w’ create a new file (default) ‘a’ append an image to a file or creates a new one if the file does
not exist (only supported by TIFF)
Supported file formats are listed by the function vigra.impexListFormats(). The different file formats support the following pixel types:
- BMP:
- Microsoft Windows bitmap image file (pixel type: UINT8 as gray and RGB).
- GIF:
- CompuServe graphics interchange format; 8-bit color (pixel type: UINT8 as gray and RGB).
- JPEG:
- Joint Photographic Experts Group JFIF format; compressed 24-bit color (pixel types: UINT8 as gray and RGB). (only available if libjpeg is installed)
- PNG:
- Portable Network Graphic (pixel types: UINT8 and UINT16 with up to 4 channels). (only available if libpng is installed)
- PBM:
- Portable bitmap format (black and white).
- PGM:
- Portable graymap format (pixel types: UINT8, INT16, INT32 as gray scale)).
- PNM:
- Portable anymap (pixel types: UINT8, INT16, INT32, gray and RGB)
- PPM:
- Portable pixmap format (pixel types: UINT8, INT16, INT32 as RGB)
- SUN:
- SUN Rasterfile (pixel types: UINT8 as gray and RGB).
- TIFF:
- Tagged Image File Format (pixel types: UINT8, INT16, INT32, FLOAT, DOUBLE with up to 4 channels). (only available if libtiff is installed.)
- VIFF:
- Khoros Visualization image file (pixel types: UINT8, INT16 INT32, FLOAT, DOUBLE with arbitrary many channels).
Wrtie a volume as a sequence of images:
writeVolume(volume, filename_base, filename_ext, dtype = '', compression = '')
The resulting image sequence will be enumerated in the form:
filename_base+[0-9]+filename_ext
Parameters ‘dtype’ and ‘compression’ will be handled as in writeImage().
vigranumpy supports all arithmetic and algebraic functions defined in numpy.ufunc, but re-implements them in module vigra.ufunc to take full advantage of axistags.
The following mathematical functions are available in this module (refer to numpy for detailed documentation):
absolute absolute add arccos arccosh arcsin arcsinh arctan
arctan2 arctanh bitwise_and bitwise_or bitwise_xor ceil
conjugate conjugate cos cosh deg2rad degrees divide equal
exp exp2 expm1 fabs floor floor_divide fmax fmin
fmod frexp greater greater_equal hypot invert invert isfinite
isinf isnan ldexp left_shift less less_equal log log10
log1p logaddexp logaddexp2 logical_and logical_not logical_or
logical_xor maximum minimum modf multiply negative not_equal
ones_like power rad2deg radians reciprocal remainder remainder
right_shift rint sign signbit sin sinh sqrt square
subtract tan tanh true_divide trunc
Some of these functions are also provided as member functions of VigraArray:
__abs__ __add__ __and__ __div__ __divmod__ __eq__
__floordiv__ __ge__ __gt__ __invert__ __le__ __lshift__
__lt__ __mod__ __mul__ __ne__ __neg__ __or__ __pos__
__pow__ __radd__ __radd__ __rand__ __rdiv__ __rdivmod__
__rfloordiv__ __rlshift__ __rmod__ __rmul__ __ror__ __rpow__
__rrshift__ __rshift__ __rsub__ __rtruediv__ __rxor__ __sub__
__truediv__ __xor__
As usual, these functions are applied independently at each pixel.
Vigranumpy overloads the numpy-versions of these functions in order to make their behavior more suitable for image analysis. In particular, we changed two aspects:
Default output types are thus determined according to the following rules:
The output type does not depend on the order of the arguments:
a + b results in the same type as b + a
- 2.a With exception of logical functions and abs(), the output type
- does not depend on the function to be executed.
2.b The output type of logical functions is bool.
- 2.c The output type of abs() follows general rules unless the
input contains complex numbers, in which case the output type is the corresponding float number type:
a + b results in the same type as a / b a == b => bool abs(complex128) => float64
If the inputs have the same type, the type is preserved:
uint8 + uint8 => uint8If (and only if) one of the inputs has at least 64 bits, the output will also have at least 64 bits:
int64 + uint32 => int64 int64 + 1.0 => float64If an array is combined with a scalar of the same kind (integer, float, or complex), the array type is preserved. If an integer array with at most 32 bits is combined with a float scalar, the result is float32 (and rule 4 kicks in if the array has 64 bits):
uint8 + 1 => uint8 uint8 + 1.0 => float32 float32 + 1.0 => float32 float64 + 1.0 => float64Integer expressions with mixed types always produce signed results. If the arguments have at most 32 bits, the result will be int32, otherwise it will be int64 (cf. rule 4):
int8 + uint8 => int32 int32 + uint8 => int32 int32 + uint32 => int32 int32 + int64 => int64 int64 + uint64 => int64In all other cases, the output type is equal to the highest input type:
int32 + float32 => float32 float32 + complex128 => complex128All defaults can be overridden by providing an explicit output array:
ufunc.add(uint8, uint8, uint16) => uint16
In order to prevent overflow, necessary upcasting is performed before the function is executed.
The module vigra.colors provides functions to adjust image brightness and contrast, and to transform between different color spaces. See Color Conversions in the C++ documentation for more information.
Adjust the brightness of a 2D scalar or multiband image. The function applies the formula:
out = image + 0.25 * log(factor) * (range[1] - range[0])
to each element of the array. ‘factor’ and ‘range[1] - range[0]’ must be positive. Elements outside the given range are clipped at the range borders. If ‘range’ is None or “” or “auto”, the range is set to the actual range of ‘image’:
range = image.min(), image.max()
Adjust the contrast of an image or volume. The function applies the formula:
out = factor * image + (1.0 - factor) * (range[1] - range[0]) / 2.0
to each element of the array. ‘factor’ and ‘range[1] - range[0]’ must be positive. Elements outside the given range are clipped at the range borders. If ‘range’ is None or “” or “auto”, the range is set to the actual range of ‘image’:
range = image.min(), image.max()
Adjust gamma correction to an image or volume. The function applies the formula:
diff = range[1] - range[0]
out = pow((image - range[0]) / diff, 1.0 / gamma) * diff + range[0]
to each element of the array. ‘gamma’ and ‘range[1] - range[0]’ must be positive. Elements outside the given range are clipped at the range borders. If ‘range’ is None or “” or “auto”, the range is set to the actual range of ‘image’:
range = image.min(), image.max()
Convert the intensity range of a 2D scalar or multiband image. The function applies a linear transformation to the intensities such that the value oldRange[0] is mapped onto newRange[0], and oldRange[1] is mapped onto newRange[1]. That is, the algorithm applies the formula:
oldDiff = oldRange[1] - oldRange[0]
newDiff = newRange[1] - newRange[0]
out = (image - oldRange[0]) / oldDiff * newDiff + newRange[0]
to each element of the array. ‘oldDiff’ and ‘newDiff’ must be positive. If ‘oldRange’ is None or “” or “auto” (the default), the range is set to the actual range of ‘image’:
range = image.min(), image.max()
If ‘newRange’ is None or “” or “auto”, it is set to (0, 255.0). If ‘out’ is explicitly passed, it must be a uin8 image.
Convert the colors of the given ‘image’ using Lab2RGBFunctor.
For details see Lab2RGBFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using Lab2RGBPrimeFunctor.
For details see Lab2RGBPrimeFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using Lab2XYZFunctor.
For details see Lab2XYZFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using Luv2RGBFunctor.
For details see Luv2RGBFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using Luv2RGBPrimeFunctor.
For details see Luv2RGBPrimeFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using Luv2XYZFunctor.
For details see Luv2XYZFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using RGB2LabFunctor.
For details see RGB2LabFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using RGB2LuvFunctor.
For details see RGB2LuvFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using RGB2RGBPrimeFunctor.
For details see RGB2RGBPrimeFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using RGB2XYZFunctor.
For details see RGB2XYZFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using RGB2sRGBFunctor.
For details see RGB2sRGBFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using RGBPrime2LabFunctor.
For details see RGBPrime2LabFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using RGBPrime2LuvFunctor.
For details see RGBPrime2LuvFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using RGBPrime2RGBFunctor.
For details see RGBPrime2RGBFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using RGBPrime2XYZFunctor.
For details see RGBPrime2XYZFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using RGBPrime2YPrimeCbCrFunctor.
For details see RGBPrime2YPrimeCbCrFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using RGBPrime2YPrimeIQFunctor.
For details see RGBPrime2YPrimeIQFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using RGBPrime2YPrimePbPrFunctor.
For details see RGBPrime2YPrimePbPrFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using RGBPrime2YPrimeUVFunctor.
For details see RGBPrime2YPrimeUVFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using XYZ2LabFunctor.
For details see XYZ2LabFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using XYZ2LuvFunctor.
For details see XYZ2LuvFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using XYZ2RGBFunctor.
For details see XYZ2RGBFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using XYZ2RGBPrimeFunctor.
For details see XYZ2RGBPrimeFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using YPrimeCbCr2RGBPrimeFunctor.
For details see YPrimeCbCr2RGBPrimeFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using YPrimeIQ2RGBPrimeFunctor.
For details see YPrimeIQ2RGBPrimeFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using YPrimePbPr2RGBPrimeFunctor.
For details see YPrimePbPr2RGBPrimeFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using YPrimeUV2RGBPrimeFunctor.
For details see YPrimeUV2RGBPrimeFunctor in the C++ documentation.
Convert the colors of the given ‘image’ using sRGB2RGBFunctor.
For details see sRGB2RGBFunctor in the C++ documentation.
The module vigra.filters provides operators that consider a window around each pixel, compute one or several numbers from the values in the window, and store the results in the corresponding pixel of the output image. This includes convolution, non-linear diffusion, morphological operators, feature detectors (such as the structure tensor) etc.
Generic 1 dimensional convolution kernel.
This kernel may be used for convolution of 1 dimensional signals or for separable convolution of multidimensional signals. The kernel’s size is given by its left() and right() methods. The desired border treatment mode is returned by getBorderTreatment(). The different init functions create a kernel with the specified properties. For more details, see Kernel1D in the C++ documentation.
Standard constructor:
Kernel1D()
Creates an identity kernel.
Copy constructor:
Kernel1D(other_kernel)
Init kernel as an averaging filter with given radius (i.e. window size 2*radius+1). ‘norm’ denotes the sum of all bins of the kernel.
Kernel construction and initialization can be performed in one step by calling the factory function ‘averagingKernel()’.
Init kernel as a binomial filter with given radius (i.e. window size 2*radius+1). ‘norm’ denotes the sum of all bins of the kernel.
Kernel construction and initialization can be performed in one step by calling the factory function ‘binomialKernel()’.
Init kernel as a 5-tap smoothing filter of the form:
[ a, 0.25, 0.5 - 2*a, 0.25, a]
Kernel construction and initialization can be performed in one step by calling the factory function ‘burtFilterKernel()’.
Init kernel as Lindeberg’s discrete analog of the Gaussian function. The radius of the kernel is always 3*std_dev. ‘norm’ denotes the desired sum of all bins of the kernel.
Kernel construction and initialization can be performed in one step by calling the factory function ‘discreteGaussianKernel()’.
Init the kernel with explicit values from ‘contents’, which must be a 1D numpy.ndarray. ‘left’ and ‘right’ are the boundaries of the kernel (inclusive). If ‘contents’ contains the wrong number of values, a run-time error results. It is, however, possible to give just one initializer. This creates an averaging filter with the given constant. The norm is set to the sum of the initializer values.
Kernel construction and initialization can be performed in one step by calling the factory function ‘explicitlyKernel()’.
Init kernel as a sampled Gaussian function. The radius of the kernel is always 3*std_dev. ‘norm’ denotes the desired sum of all bins of the kernel (i.e. the kernel is corrected for the normalization error introduced by windowing the Gaussian to a finite interval). However, if norm is 0.0, the kernel is normalized to 1 by the analytic expression for the Gaussian, and no correction for the windowing error is performed.
Kernel construction and initialization can be performed in one step by calling the factory function ‘gaussianKernel()’.
Init kernel as a Gaussian derivative of order ‘order’. The radius of the kernel is always 3*std_dev + 0.5*order. ‘norm’ denotes the norm of the kernel. Thus, the kernel will be corrected for the error introduced by windowing the Gaussian to a finite interval. However, if norm is 0.0, the kernel is normalized to 1 by the analytic expression for the Gaussian derivative, and no correction for the windowing error is performed.
Kernel construction and initialization can be performed in one step by calling the factory function ‘gaussianDerivativeKernel()’.
Init kernel as a 3-tap second difference filter of the form:
[ 1, -2, 1]
Kernel construction and initialization can be performed in one step by calling the factory function ‘secondDifference3Kernel()’.
Init kernel as a symmetric difference filter of the form:
[ 0.5 * norm, 0.0 * norm, -0.5 * norm]
Kernel construction and initialization can be performed in one step by calling the factory function ‘symmetricDifferenceKernel()’.
Generic 2 dimensional convolution kernel.
This kernel may be used for convolution of 2 dimensional signals. The desired border treatment mode is returned by borderTreatment().(Note that the 2D convolution functions don’t currently support all modes.) The different init functions create a kernel with the specified properties. For more details, see Kernel2D in the C++ documentation.
Standard constructor:
Kernel2D()
Creates an identity kernel.
Copy constructor:
Kernel2D(other_kernel)
Init the 2D kernel as a circular averaging filter. The norm will be calculated as 1 / (number of non-zero kernel values).
Precondition:
radius > 0
Kernel construction and initialization can be performed in one step by calling the factory function ‘diskKernel2D()’.
Init the kernel with explicit values from ‘contents’, which must be a 2D numpy.ndarray. ‘upperLeft’ and ‘lowerRight’ are the boundaries of the kernel (inclusive), and must be 2D tuples. If ‘contents’ contains the wrong number of values, a run-time error results. It is, however, possible to give just one initializer. This creates an averaging filter with the given constant. The norm is set to the sum of the initializer values.
Kernel construction and initialization can be performed in one step by calling the factory function ‘explicitlyKernel2D()’.
Init kernel as a sampled 2D Gaussian function. The radius of the kernel is always 3*std_dev. ‘norm’ denotes the desired sum of all bins of the kernel (i.e. the kernel is corrected for the normalization error introduced by windowing the Gaussian to a finite interval). However, if norm is 0.0, the kernel is normalized to 1 by the analytic expression for the Gaussian, and no correction for the windowing error is performed.
Kernel construction and initialization can be performed in one step by calling the factory function ‘gaussianKernel2D()’.
Init the 2D kernel as the cartesian product of two 1D kernels of type Kernel1D. The norm becomes the product of the two original norms.
Kernel construction and initialization can be performed in one step by calling the factory function ‘separableKernel2D()’.
Convolve an image with the given ‘kernel’ (or kernels). If the input has multiple channels, the filter is applied to each channel independently. The function can be used in 3 different ways:
For details see separableConvolveMultiArray and convolveImage in the vigra C++ documentation.
Convolution along a single dimension of a 2D scalar or multiband image. ‘kernel’ must be an instance of Kernel1D.
For details see convolveMultiArrayOneDimension in the vigra C++ documentation.
Apply a closing filter with disc of given radius to image.
This is an abbreviation for applying a dilation and an erosion filter in sequence. This function also works for multiband images, it is then executed on every band.
See discRankOrderFilter in the C++ documentation for more information.
Apply dilation (maximum) filter with disc of given radius to image.
This is an abbreviation for the rank order filter with rank = 1.0. This function also works for multiband images, it is then executed on every band.
See discDilation in the C++ documentation for more information.
Apply erosion (minimum) filter with disc of given radius to image.
This is an abbreviation for the rank order filter with rank = 0.0. This function also works for multiband images, it is then executed on every band.
See discErosion in the C++ documentation for more information.
Apply median filter with disc of given radius to image.
This is an abbreviation for the rank order filter with rank = 0.5. This function also works for multiband images, it is then executed on every band.
See discMedian in the C++ documentation for more information.
Apply a opening filter with disc of given radius to image.
This is an abbreviation for applying an erosion and a dilation filter in sequence. This function also works for multiband images, it is then executed on every band.
See discRankOrderFilter in the C++ documentation for more information.
Apply rank order filter with disc structuring function to a float image.
The pixel values of the source image must be in the range 0...255. Radius must be >= 0. Rank must be in the range 0.0 <= rank <= 1.0. The filter acts as a minimum filter if rank = 0.0, as a median if rank = 0.5, and as a maximum filter if rank = 1.0. This function also works for multiband images, it is then executed on every band.
For details see discRankOrderFilter in the C++ documentation.
Apply rank order filter with disc structuring function to a float image using a mask.
The pixel values of the source image must be in the range 0...255. Radius must be >= 0.Rank must be in the range 0.0 <= rank <= 1.0. The filter acts as a minimum filter if rank = 0.0,as a median if rank = 0.5, and as a maximum filter if rank = 1.0.
The mask is only applied to the input image, i.e. the function generates an output wherever the current disc contains at least one pixel with non-zero mask value. Source pixels with mask value zero are ignored during the calculation of the rank order.
This function also works for multiband images, it is then executed on every band. If the mask has only one band, it is used for every image band. If the mask has the same number of bands, as the image the bands are used for the corresponding image bands.
For details see discRankOrderFilterWithMask in the C++ documentation.
Compute the distance transform of a 2D scalar float image. All pixels with a value of 0.0 are considered to be background pixels, while all pixels with a nonzero value are considered to be foreground pixels. The parameter ‘background’ is a Boolean scalar that specifies whether to compute the distance of all background pixels to the nearest foreground pixels (if it is ‘True’, default) or vice versa (if it is ‘False’). Hence in the destination image, for background==True all background pixels will be assigned their distance value, while all foreground pixels will be assigned 0. For background==False, it is exactly the other way around.
The ‘norm’ parameter gives the distance norm to use (0: infinity norm, 1: L1 norm, 2: Euclidean norm).
For details see distanceTransform in the vigra C++ documentation.
Compute the Euclidean distance transform of a 3D scalar float volume. All voxels with a value of 0.0 are considered to be background voxels, while all voxels with a nonzero value are considered to be foreground voxels. The parameter ‘background’ is a Boolean scalar that specifies whether to compute the distance of all background voxels to the nearest foreground voxel (if it is ‘True’, default) or vice versa (if it is ‘False’). Hence in the destination volume, for background==True all background voxels will be assigned their distance value, while all foreground voxels will be assigned 0. For background==False, it is exactly the other way around.
For more details see separableMultiDistance in the vigra C++ documentation.
Calculate the gradient vector by means of a 1st derivative of Gaussian filter at the given scale for a 2D scalar image.
If ‘sigma’ is a single value, an isotropic filter at this scale is applied (i.e., each dimension is filtered in the same way). If ‘sigma’ is a tuple or list of values, the amount of smoothing will be different for each spatial dimension. The optional ‘sigma_d’ (single, tuple, or list) denotes the resolution standard deviation per axis, the optional ‘step_size’ (single, tuple, or list) the distance between two adjacent pixels for each dimension. The length of the tuples or lists must be equal to the number of spatial dimensions.
‘window_size’ and ‘roi’ have the same meaning as in gaussianSmoothing().
For details see gaussianGradientMultiArray and ConvolutionOptions in the vigra C++ documentation.
Calculate the gradient magnitude by means of a 1st derivative of Gaussian filter at the given scale for a 2D scalar or multiband image. If ‘accumulate’ is True (the default), the gradients are accumulated (in the L2-norm sense) over all channels of a multi-channel array. Otherwise, a separate gradient magnitude is computed for each channel.
If ‘sigma’ is a single value, an isotropic filter at this scale is applied (i.e., each dimension is filtered in the same way). If ‘sigma’ is a tuple or list of values, the amount of smoothing will be different for each spatial dimension. The optional ‘sigma_d’ (single, tuple, or list) denotes the resolution standard deviation per axis, the optional ‘step_size’ (single, tuple, or list) the distance between two adjacent pixels for each dimension. The length of the tuples or lists must be equal to the number of spatial dimensions.
‘window_size’ and ‘roi’ have the same meaning as in gaussianSmoothing().
For details see gaussianGradientMultiArray and ConvolutionOptions in the vigra C++ documentation.
Perform sharpening function with gaussian filter.
For details see gaussianSharpening in the vigra C++ documentation.
Perform Gaussian smoothing of a 2D or 3D scalar or multiband array.
Each channel of the array is smoothed independently. If ‘sigma’ is a single value, an isotropic Gaussian filter at this scale is applied (i.e. each dimension is smoothed in the same way). If ‘sigma’ is a tuple or list of values, the amount of smoothing will be different for each spatial dimension. The optional ‘sigma_d’ (single, tuple, or list) (single, tuple, or list) denotes the resolution standard deviation per axis, the optional ‘step_size’ (single, tuple, or list) the distance between two adjacent pixels for each dimension. The length of the tuples or lists must be equal to the number of spatial dimensions.
‘window_size’ specifies the ratio between the effective filter scale and the size of the filter window. Use a value around 2.0 to speed-up the computation by increasing the error resulting from cutting off the Gaussian. For the default 0.0, the window size is automatically determined.
If ‘roi’ is not None, it must specify the desired region-of-interest as a pair ‘(first_point, beyond_last_point)’ (e.g. ‘roi=((10,20), (200,250))’). As usual, the second point is the first point outside the ROI, and the ROI must not be outside the input array dimensions. The coordinates refer only to non-channel axes - if your array has an explicit channel axis, the ROI dimension must be one less than the array dimension. If you pass in an explicit ‘out’ array and specify an ROI, the ‘out’ array must have the shape of the ROI.
For details see gaussianSmoothing and ConvolutionOptions in the vigra C++ documentation.
Calculate the Hessian matrix by means of a derivative of Gaussian filters at the given scale for a 2D scalar image.
If ‘sigma’ is a single value, an isotropic filter at this scale is applied (i.e., each dimension is filtered in the same way). If ‘sigma’ is a tuple or list of values, the amount of smoothing will be different for each spatial dimension. The optional ‘sigma_d’ (single, tuple, or list) denotes the resolution standard deviation per axis, the optional ‘step_size’ (single, tuple, or list) the distance between two adjacent pixels for each dimension. The length of the tuples or lists must be equal to the number of spatial dimensions.
‘window_size’ and ‘roi’ have the same meaning as in gaussianSmoothing().
For details see hessianOfGaussianMultiArray in the vigra C++ documentation.
Calculate the Hessian matrix by means of a derivative of Gaussian filters at the given scale for a 2D scalar image.
If ‘sigma’ is a single value, an isotropic filter at this scale is applied (i.e., each dimension is filtered in the same way). If ‘sigma’ is a tuple or list of values, the amount of smoothing will be different for each spatial dimension. The optional ‘sigma_d’ (single, tuple, or list) denotes the resolution standard deviation per axis, the optional ‘step_size’ (single, tuple, or list) the distance between two adjacent pixels for each dimension. The length of the tuples or lists must be equal to the number of spatial dimensions.
‘window_size’ and ‘roi’ have the same meaning as in gaussianSmoothing().
For details see hessianOfGaussianMultiArray and ConvolutionOptions in the vigra C++ documentation.
Calculate the Hessian matrix by means of a derivative of Gaussian filters at the given scale for a 3D scalar image.
For details see hessianOfGaussianMultiArray in the vigra C++ documentation.
Compute the eigenvalues of the Hessian of Gaussian at the given scale for a scalar image or volume.
Calls hessianOfGaussian() and tensorEigenvalues().
Anisotropic tensor smoothing with the hourglass filter.
For details see hourGlassFilter in the vigra C++ documentation.
Filter 2D or 3D scalar array with the Laplacian of Gaussian operator at the given scale.
If ‘sigma’ is a single value, an isotropic filter at this scale is applied (i.e., each dimension is filtered in the same way). If ‘sigma’ is a tuple or list of values, the amount of smoothing will be different for each spatial dimension. The optional ‘sigma_d’ (single, tuple, or list) denotes the resolution standard deviation per axis, the optional ‘step_size’ (single, tuple, or list) the distance between two adjacent pixels for each dimension. The length of the tuples or lists must be equal to the number of spatial dimensions.
‘window_size’ and ‘roi’ have the same meaning as in gaussianSmoothing().
For details see laplacianOfGaussianMultiArray and ConvolutionOptions in the vigra C++ documentation.
Binary closing on a 3D scalar or multiband uint8 array.
This function applies a flat circular opening operator (sequential dilation and erosion) with a given radius. The operation is isotropic. The input is a uint8 or boolean multi-dimensional array where non-zero pixels represent foreground and zero pixels represent background. This function also works for multiband arrays, it is then executed on every band.
For details see vigra C++ documentation (multiBinaryDilation and multiBinaryErosion).
Binary dilation on a 3D scalar or multiband uint8 array.
This function applies a flat circular dilation operator with a given radius. The operation is isotropic. The input is a uint8 or boolean multi-dimensional array where non-zero pixels represent foreground and zero pixels represent background. This function also works for multiband arrays, it is then executed on every band.
For details see multiBinaryDilation in the C++ documentation.
Binary erosion on a 3D scalar or multiband uint8 array.
This function applies a flat circular erosion operator with a given radius. The operation is isotropic. The input is a uint8 or boolean multi-dimensional array where non-zero pixels represent foreground and zero pixels represent background. This function also works for multiband arrays, it is then executed on every band.
For details see multiBinaryErosion in the C++ documentation.
Binary opening on a 3D scalar or multiband uint8 array.
This function applies a flat circular opening operator (sequential erosion and dilation) with a given radius. The operation is isotropic. The input is a uint8 or boolean multi-dimensional array where non-zero pixels represent foreground and zero pixels represent background. This function also works for multiband arrays, it is then executed on every band.
For details see vigra C++ documentation (multiBinaryDilation and multiBinaryErosion).
Parabolic grayscale closing on multi-dimensional arrays.
This function applies a parabolic closing (sequential dilation and erosion) operator with a given spread ‘sigma’ on a grayscale array. The operation is isotropic. The input is a grayscale multi-dimensional array. This function also works for multiband arrays, it is then executed on every band.
For details see multiGrayscaleDilation and multiGrayscaleErosion in the C++ documentation.
multiGrayscaleClosing( (object)image, (float)sigma [, (object)out=None]) -> object
Parabolic grayscale dilation on multi-dimensional arrays.
This function applies a parabolic dilation operator with a given spread ‘sigma’ on a grayscale array. The operation is isotropic. The input is a grayscale multi-dimensional array. This function also works for multiband arrays, it is then executed on every band.
For details see multiGrayscaleDilation in the C++ documentation.
multiGrayscaleDilation( (object)image, (float)sigma [, (object)out=None]) -> object
Parabolic grayscale erosion on a 3D scalar or multiband uint8 array.
This function applies a parabolic erosion operator with a given spread ‘sigma’ on a grayscale array. The operation is isotropic. The input is a grayscale multi-dimensional array. This function also works for multiband arrays, it is then executed on every band.
For details see multiGrayscaleErosion in the C++ documentation.
Parabolic grayscale opening on multi-dimensional arrays.
This function applies a parabolic opening (sequential erosion and dilation) operator with a given spread ‘sigma’ on a grayscale array. The operation is isotropic. The input is a grayscale multi-dimensional array. This function also works for multiband arrays, it is then executed on every band.
For details see multiGrayscaleDilation and multiGrayscaleErosion in the C++ documentation.
multiGrayscaleOpening( (object)image, (float)sigma [, (object)out=None]) -> object
Perform edge-preserving smoothing at the given scale.
For details see nonlinearDiffusion in the vigra C++ documentation.
Perform normalized convolution of an image. If the image has multiple channels, every channel is convolved independently. The ‘mask’ tells the algorithm whether input pixels are valid (non-zero mask value) or not. Invalid pixels are ignored in the convolution. The mask must have one channel (which is then used for all channels input channels) or as many channels as the input image.
For details, see normalizedConvolveImage in the C++ documentation.
Find centers of radial symmetry in an 2D image.
This algorithm implements the Fast Radial Symmetry Transform according to [G. Loy, A. Zelinsky: “A Fast Radial Symmetry Transform for Detecting Points of Interest”, in: A. Heyden et al. (Eds.): Proc. of 7th European Conf. on Computer Vision, Part 1, pp. 358-368, Springer LNCS 2350, 2002]
For details see radialSymmetryTransform in the vigra C++ documentation.
Perform 2D convolution with a first-order recursive filter with parameter ‘b’ and given ‘borderTreatment’. ‘b’ must be between -1 and 1.
For details see recursiveFilterX and recursiveFilterY (which this function calls in succession) in the vigra C++ documentation.
Perform 2D convolution with a second-order recursive filter with parameters ‘b1’ and ‘b2’. Border treatment is always BORDER_TREATMENT_REFLECT.
For details see recursiveFilterX and recursiveFilterY (which this function calls in succession) in the vigra C++ documentation.
Compute a fast approximate Gaussian smoothing of a 2D scalar or multiband image.
This function uses the third-order recursive filter approximation to the Gaussian filter proposed by Young and van Vliet. Each channel of the array is smoothed independently. If ‘sigma’ is a single value, an isotropic Gaussian filter at this scale is applied (i.e. each dimension is smoothed in the same way). If ‘sigma’ is a tuple of values, the amount of smoothing will be different for each spatial dimension. The length of the tuple must be equal to the number of spatial dimensions.
For details see recursiveGaussianFilterLine in the vigra C++ documentation.
Compute the gradient of a scalar image using a recursive (exponential) filter at the given ‘scale’. The output image (if given) must have two channels.
For details see recursiveSmoothLine and recursiveFirstDerivativeLine (which this function calls internally) in the vigra C++ documentation.
Compute the gradient of a 2D scalar or multiband image using a recursive (exponential) filter at the given ‘scale’. The output image (if given) must have as many channels as the input.
For details see recursiveSmoothLine and recursiveSecondDerivativeLine (which this function calls internally) in the vigra C++ documentation.
Calls recursiveFilter2D() with b = exp(-1/scale), which corresponds to smoothing with an exponential filter exp(-abs(x)/scale).
For details see recursiveSmoothLine in the vigra C++ documentation.
Calculate Riesz transforms of the Laplacian of Gaussian.
For details see rieszTransformOfLOG in the vigra C++ documentation.
Perform simple sharpening function.
For details see simpleSharpening in the vigra C++ documentation.
Calculate the structure tensor of an image by means of Gaussian (derivative) filters at the given scales. If the input has multiple channels, the structure tensors of each channel are added to get the result.
If ‘innerScale’ and ‘outerScale’ are single values, isotropic filters at these scales are applied (i.e., each dimension is filtered in the same way). If ‘innerScale’ and / or ‘outerScale’ are are tuples or lists of values, the amount of smoothing will be different for each spatial dimension. The optional ‘sigma_d’ (single, tuple, or list) denotes the resolution standard deviation per axis, the optional ‘step_size’ (single, tuple, or list) the distance between two adjacent pixels for each dimension. The length of the tuples or lists must be equal to the number of spatial dimensions.
‘window_size’ and ‘roi’ have the same meaning as in gaussianSmoothing().
For details see structureTensorMultiArray and ConvolutionOptions in the vigra C++ documentation.
Compute the eigenvalues of the structure tensor at the given scales for a scalar or multi-channel image or volume.
Calls structureTensor() and tensorEigenvalues().
Calculate gradient of a scalar 2D image using symmetric difference filters. The optional tuple or list ‘step_size’ denotes the distance between two adjacent pixels for each dimension; its length must be equal to the number of spatial dimensions.
‘roi’ has the same meaning as in gaussianSmoothing().
For details see symmetricGradientMultiArray and ConvolutionOptions in the vigra C++ documentation.
Calculate the determinant of a 2x2 tensor image.
For details see tensorDeterminantMultiArray in the vigra C++ documentation.
Calculate eigen representation of a symmetric 2x2 tensor.
For details see tensorEigenRepresentation in the vigra C++ documentation.
Calculate the eigenvalues in each pixel/voxel of a 2x2 tensor image.
For details see tensorEigenvaluesMultiArray in the vigra C++ documentation.
Calculate the trace of a 2x2 tensor image.
For details see tensorTraceMultiArray in the vigra C++ documentation.
Turn a 2D vector valued image (e.g. the gradient image) into a tensor image by computing the outer product in every pixel.
For details see vectorToTensorMultiArray in the vigra C++ documentation.
The module vigra.sampling contains methods to change the number and/or location of the image sampling points, such as resizing, rotation, and interpolation.
Resample an image by the given ‘factor’
The ‘out’ parameter must have, if given, the according dimensions. This function also works for multiband images, it is then executed on every band.
For more details, see resampleImage in the vigra C++ documentation.
Resample image using a gaussian filter:
resamplingGaussian(image,
sigmaX=1.0, derivativeOrderX=0, samplingRatioX=2.0, offsetX=0.0,
sigmaY=1.0, derivativeOrderY=0, samplingRatioY=2.0, offsetY=0.0,
out=None)
This function utilizes resamplingConvolveImage with a Gaussianfilter (see the vigra C++ documentation for details).
Resize image or volume using B-spline interpolation.
The spline order is given in the parameter ‘order’. The desired shape of the output array is taken either from ‘shape’ or ‘out’. If both are given, they must agree. This function also works for multi-channel data, it is then executed on every channel independently.
For more details, see resizeImageSplineInterpolation and resizeMultiArraySplineInterpolation in the vigra C++ documentation.
resize( (object)image [, (object)shape=None [, (int)order=3 [, (object)out=None]]]) -> object
Resize image using the Catmull/Rom interpolation function.
The desired shape of the output image is taken either from ‘shape’ or ‘out’. If both are given, they must agree. This function also works for multiband images, it is then executed on every band.
For more details, see resizeImageCatmullRomInterpolation in the vigra C++ documentation.
Resize image using the Coscot interpolation function.
The desired shape of the output image is taken either from ‘shape’ or ‘out’. If both are given, they must agree. This function also works for multiband images, it is then executed on every band.
For more details, see resizeImageCoscotInterpolation in the vigra C++ documentation.
Resize image using linear interpolation. The function uses the standard separable bilinear interpolation algorithm to obtain a good compromise between quality and speed.
The desired shape of the output image is taken either from ‘shape’ or ‘out’. If both are given, they must agree. This function also works for multiband images, it is then executed on every band.
For more details, see resizeImageLinearInterpolation in the vigra C++ documentation.
Resize image by repeating the nearest pixel values.
The desired shape of the output image is taken either from ‘shape’ or ‘out’. If both are given, they must agree. This function also works for multiband images, it is then executed on every band.
For more details, see resizeImageNoInterpolation in the vigra C++ documentation.
Resize image using B-spline interpolation.
The spline order is given in the parameter ‘order’. The desired shape of the output image is taken either from ‘shape’ or ‘out’. If both are given, they must agree. This function also works for multiband images, it is then executed on every band.
For more details, see resizeImageSplineInterpolation in the vigra C++ documentation.
Resize volume using B-spline interpolation.
The spline order is given in the parameter ‘order’. The dimensions of the output volume is taken either from ‘shape’ or ‘out’. If both are given, they must agree. This function also works for multiband volumes, it is then executed on every band.
For more details, see resizeMultiArraySplineInterpolation in the vigra C++ documentation.
Rotate an image by an arbitrary angle using splines for interpolation around its center.
The angle may be given in degree (parameter degree). The parameter ‘splineOrder’ indicates the order of the splines used for interpolation. If the ‘out’ parameter is given, the image is cropped for it’s dimensions. If the ‘out’ parameter is not given, an output image with the same dimensions as the input image is created.
For more details, see GeometricTransformations.rotationMatrix2DDegrees in the vigra C++ documentation.
Rotate an image by an arbitrary angle around its center using splines for interpolation.
The angle may be given in radiant (parameter radiant). The parameter ‘splineOrder’ indicates the order of the splines used for interpolation. If the ‘out’ parameter is given, the image is cropped for it’s dimensions. If the ‘out’ parameter is not given, an output image with the same dimensions as the input image is created.
For more details, see GeometricTransformations.rotationMatrix2DRadians in the vigra C++ documentation.
Rotate an image by a multiple of 90 degrees.
The ‘orientation’ parameter (which must be one of CLOCKWISE, COUNTER_CLOCKWISE and UPSIDE_DOWN indicates the rotation direction. The ‘out’ parameter must, if given, have the according dimensions. This function also works for multiband images, it is then executed on every band.
For more details, see rotateImage in the vigra C++ documentation.
Spline image views implement an interpolated view for an image which can be accessed at real-valued coordinates (in contrast to the plain image, which can only be accessed at integer coordinates). Module vigra.sampling defines:
SplineImageView0
SplineImageView1
SplineImageView2
SplineImageView3
SplineImageView4
SplineImageView5
The number denotes the spline interpolation order of the respective classes. Below, we describe SplineImageView3 in detail, but the other classes work analogously. See SplineImageView in the C++ documentation for more detailed information.
Construct a SplineImageView for the given image:
SplineImageView(image, skipPrefilter = False)
Currently, ‘image’ can have dtype numpy.uint8, numpy.int32, and numpy.float32. If ‘skipPrefilter’ is True, image values are directly used as spline coefficients, so that the view performs approximation rather than interploation.
__init__( (object)arg1, (object)arg2) -> object
__init__( (object)arg1, (object)arg2) -> object
__init__( (object)arg1, (object)arg2, (bool)arg3) -> object
__init__( (object)arg1, (object)arg2, (bool)arg3) -> object
__init__( (object)arg1, (object)arg2, (bool)arg3) -> object
Return first derivative in x direction at a real-valued coordinate.
SplineImageView.dx(x, y) -> value
Return third derivative in x direction at a real-valued coordinate.
SplineImageView.dx3(x, y) -> value
Like dx3(), but returns an entire image with the given sampling factors. For example,
SplineImageView.dx3Image(2.0, 2.0) -> image
creates an derivative image with two-fold oversampling in both directions.
Like dx(), but returns an entire image with the given sampling factors. For example,
SplineImageView.dxImage(2.0, 2.0) -> image
creates an derivative image with two-fold oversampling in both directions.
Return second derivative in x direction at a real-valued coordinate.
SplineImageView.dxx(x, y) -> value
Like dxx(), but returns an entire image with the given sampling factors. For example,
SplineImageView.dxxImage(2.0, 2.0) -> image
creates an derivative image with two-fold oversampling in both directions.
Return mixed third derivative at a real-valued coordinate.
SplineImageView.dxxy(x, y) -> value
Like dxxy(), but returns an entire image with the given sampling factors. For example,
SplineImageView.dxxyImage(2.0, 2.0) -> image
creates an derivative image with two-fold oversampling in both directions.
Return mixed second derivative at a real-valued coordinate.
SplineImageView.dxy(x, y) -> value
Like dxy(), but returns an entire image with the given sampling factors. For example,
SplineImageView.dxyImage(2.0, 2.0) -> image
creates an derivative image with two-fold oversampling in both directions.
Return mixed third derivative at a real-valued coordinate.
SplineImageView.dxyy(x, y) -> value
Like dxyy(), but returns an entire image with the given sampling factors. For example,
SplineImageView.dxyyImage(2.0, 2.0) -> image
creates an derivative image with two-fold oversampling in both directions.
Return first derivative in y direction at a real-valued coordinate.
SplineImageView.dy(x, y) -> value
Return third derivative in y direction at a real-valued coordinate.
SplineImageView.dy3(x, y) -> value
Like dy3(), but returns an entire image with the given sampling factors. For example,
SplineImageView.dy3Image(2.0, 2.0) -> image
creates an derivative image with two-fold oversampling in both directions.
Like dy(), but returns an entire image with the given sampling factors. For example,
SplineImageView.dyImage(2.0, 2.0) -> image
creates an derivative image with two-fold oversampling in both directions.
Return second derivative in y direction at a real-valued coordinate.
SplineImageView.dyy(x, y) -> value
Like dyy(), but returns an entire image with the given sampling factors. For example,
SplineImageView.dyyImage(2.0, 2.0) -> image
creates an derivative image with two-fold oversampling in both directions.
SplineImageView.facetCoefficients(x, y) -> matrix
Return the facet coefficient matrix so that spline values can be computed explicitly. The matrix has size (order+1)x(order+1), where order is the order of the spline. The matrix must be multiplied from left and right with the powers of the local facet x- and y-coordinates respectively (note that local facet coordinates are in the range [0,1] for odd order splines and [-0.5, 0.5] for even order splines).
Usage for odd spline order:
s = SplineImageView3(image) c = s.coefficients(10.1, 10.7) x = matrix([1, 0.1, 0.1**2, 0.1**3]) y = matrix([1, 0.7, 0.7**2, 0.7**3]) assert abs(x * c * y.T - s[10.1, 10.7]) < smallNumber
Usage for even spline order:
s = SplineImageView2(image) c = s.coefficients(10.1, 10.7) x = matrix([1, 0.1, 0.1**2]) y = matrix([1, -0.3, (-0.3)**2]) assert abs(x * c * y.T - s[10.1, 10.7]) < smallNumber
Return gradient squared magnitude at a real-valued coordinate.
SplineImageView.g2(x, y) -> value
Like g2(), but returns an entire image with the given sampling factors. For example,
SplineImageView.g2Image(2.0, 2.0) -> image
creates an derivative image with two-fold oversampling in both directions.
Return first derivative in x direction of the gradient squared magnitude at a real-valued coordinate.
SplineImageView.g2x(x, y) -> value
Like g2x(), but returns an entire image with the given sampling factors. For example,
SplineImageView.g2xImage(2.0, 2.0) -> image
creates an derivative image with two-fold oversampling in both directions.
Return first derivative in y direction of the gradient squared magnitude at a real-valued coordinate.
SplineImageView.g2y(x, y) -> value
Like g2y(), but returns an entire image with the given sampling factors. For example,
SplineImageView.g2yImage(2.0, 2.0) -> image
creates an derivative image with two-fold oversampling in both directions.
Return an interpolated image or derivative image with the given sampling factors and derivative orders. For example, we get a two-fold oversampled image with the x-derivatives in each pixel by:
SplineImageView.interpolatedImage(2.0, 2.0, 1, 0) -> image
Check if a coordinate is inside the underlying image.
SplineImageView.isInside(x, y) -> bool
Check if a coordinate is within the valid range of the SplineImageView.
SplineImageView.isValid(x, y) -> bool
Thanks to reflective boundary conditions, the valid range is three times as big as the size of the underlying image.
Bases: list
Create a new pyramid. The new pyramid levels range from ‘lowestLevel’ to ‘highestLevel’ (inclusive), and the given ‘image’ is copied to ‘copyImageToLevel’. The images at other levels are filled with zeros and sized so that the shape is reduced by half when going up (to higher levels), and doubled when going down.
This class can handle multi-channel images, but only when image.channelIndex exists and returns image.ndim-1 (i.e. the image must have axistags, and the channel axis must correspond to the last index, as in C- or V-order).
Expand the image at ‘srcLevel’ to ‘destLevel’, using the Burt smoothing filter with the given ‘centerValue’. srcLevel must be larger than destLevel.
For more details, see pyramidExpandBurtFilter in the C++ documentation.
Expand the image at ‘srcLevel’ to ‘destLevel’, using the Burt smoothing filter with the given ‘centerValue’, and reconstruct the images for the levels srcLevel-1 ... destLevel from their Laplacian images. srcLevel must be larger than destLevel.
For more details, see pyramidExpandBurtLaplacian in the C++ documentation.
Reduce the image at ‘srcLevel’ to ‘destLevel’, using the Burt smoothing filter with the given ‘centerValue’. srcLevel must be smaller than destLevel.
For more details, see pyramidReduceBurtFilter in the C++ documentation.
Reduce the image at ‘srcLevel’ to ‘destLevel’, using the Burt smoothing filter with the given ‘centerValue’, and compute Laplacian images for the levels srcLevel ... destLevel-1. srcLevel must be smaller than destLevel.
For more details, see pyramidReduceBurtLaplacian in the C++ documentation.
The module vigra.fourier contains functions for Fourier transforms, Cosine/Sine transforms, and Fourier-domain filters.
The module vigra.analysis contains segmentation algorithms (e.g. watershed), edge and corner detection, localization of maxima and minima etc.
Represent an Edgel at a particular subpixel position (x, y), having given ‘strength’ and ‘orientation’.
For details, see Edgel in the vigra C++ documentation.
Standard constructor:
Edgel()
Constructor:
Edgel(x, y, strength, orientation)
Beautify crack edge image for visualization.
For details see beautifyCrackEdgeImage in the vigra C++ documentation.
Detect and mark edges in an edge image using Canny’s algorithm.
For details see cannyEdgeImage in the vigra C++ documentation.
Detect and mark edges in an edge image using Canny’s algorithm.
For details see cannyEdgeImageWithThinning in the vigra C++ documentation.
Return a list of Edgel objects whose strength is at least ‘threshold’.
The function comes in two forms:
cannyEdgelList(gradient, threshold) -> list
cannyEdgelList(image, scale, threshold) -> list
The first form expects a gradient image (i.e. with two channels) to compute edgels, whereas the second form expects a scalar image and computes the gradient internally at ‘scale’.
For details see cannyEdgelList in the vigra C++ documentation.
Return a list of Edgel objects whose strength is at least ‘threshold’.
The function comes in two forms:
cannyEdgelList3x3(gradient, threshold) -> list
cannyEdgelList3x3(image, scale, threshold) -> list
The first form expects a gradient image (i.e. with two channels) to compute edgels, whereas the second form expects a scalar image and computes the gradient internally at ‘scale’. The results are slightly better than those of cannyEdgelList().
For details see cannyEdgelList3x3 in the vigra C++ documentation.
Close one-pixel wide gaps in a cell grid edge image.
For details see closeGapsInCrackEdgeImage in the vigra C++ documentation.
Find corners in a scalar 2D image using the method of Beaudet at the given ‘scale’.
For details see beaudetCornerDetector in the vigra C++ documentation.
Find corners in a scalar 2D image using the boundary tensor at the given ‘scale’.
Specifically, the cornerness is defined as twice the small eigenvalue of the boundary tensor.
For details see boundaryTensor in the vigra C++ documentation.
Find corners in a scalar 2D image using the method of Foerstner at the given ‘scale’.
For details see foerstnerCornerDetector in the vigra C++ documentation.
Find corners in a scalar 2D image using the method of Harris at the given ‘scale’.
For details see cornerResponseFunction in the vigra C++ documentation.
Find corners in a scalar 2D image using the method of Rohr at the given ‘scale’.
For details see rohrCornerDetector in the vigra C++ documentation.
Find local maxima and maximal plateaus in an image and mark them with the given ‘marker’. Parameter ‘neighborhood’ specifies the pixel neighborhood to be used and can be 4 or 8 (default).
For details see extendedLocalMaxima in the vigra C++ documentation.
Find local maxima and maximal plateaus in a volume and mark them with the given ‘marker’. Parameter ‘neighborhood’ specifies the pixel neighborhood to be used and can be 6 (default) or 26 .
For details see extendedLocalMaxima3D in the vigra C++ documentation.
Find local minima and minimal plateaus in an image and mark them with the given ‘marker’. Parameter ‘neighborhood’ specifies the pixel neighborhood to be used and can be 4 or 8 (default).
For details see extendedLocalMinima in the vigra C++ documentation.
Find local minima and minimal plateaus in a volume and mark them with the given ‘marker’. Parameter ‘neighborhood’ specifies the pixel neighborhood to be used and can be 6(default) or 26 .
For details see extendedLocalMinima3D in the vigra C++ documentation.
labelImage( (object)image [, (int)neighborhood=4 [, (object)out=None]]) -> object
Find the connected components of a segmented image. Parameter ‘neighborhood’ specifies the pixel neighborhood to be used and can be 4 (default) or 8.
For details see labelImage in the vigra C++ documentation.
labelImageWithBackground( (object)image [, (int)neighborhood=4 [, (int)background_value=0 [, (object)out=None]]]) -> object
Find the connected components of a segmented image, excluding the background from labeling, where the background is the set of all pixels with the given ‘background_value’. Parameter ‘neighborhood’ specifies the pixel neighborhood to be used and can be 4 (default) or 8.
For details see labelImageWithBackground in the vigra C++ documentation.
labelVolume( (object)volume [, (int)neighborhood=6 [, (object)out=None]]) -> object
Find the connected components of a segmented volume. Parameter ‘neighborhood’ specifies the pixel neighborhood to be used and can be 6 (default) or 26.
For details see labelVolume in the vigra C++ documentation.
labelVolumeWithBackground( (object)volume [, (int)neighborhood=6 [, (int)background_value=0 [, (object)out=None]]]) -> object
Find the connected components of a segmented volume, excluding the background from labeling, where the background is the set of all pixels with the given ‘background_value’. Parameter ‘neighborhood’ specifies the pixel neighborhood to be used and can be 6 (default) or 26.
For details see labelVolumeWithBackground in the vigra C++ documentation.
Find local maxima in an image and mark them with the given ‘marker’. Parameter ‘neighborhood’ specifies the pixel neighborhood to be used and can be 4 or 8 (default).
For details see localMaxima in the vigra C++ documentation.
Find local maxima and maximal plateaus in a volume and mark them with the given ‘marker’. Parameter ‘neighborhood’ specifies the pixel neighborhood to be used and can be 6(default) or 26 .
For details see localMaxima3D in the vigra C++ documentation.
Find local minima in an image and mark them with the given ‘marker’. Parameter ‘neighborhood’ specifies the pixel neighborhood to be used and can be 4 or 8 (default).
For details see localMinima in the vigra C++ documentation.
Find local minima in a volume and mark them with the given ‘marker’. Parameter ‘neighborhood’ specifies the pixel neighborhood to be used and can be 6 or 26 (default).
For details see localMinima3D in the vigra C++ documentation.
Transform a labeled uint32 image into a crack edge image.
For details see regionImageToCrackEdgeImage in the vigra C++ documentation.
Transform a labeled uint32 image into an edge image.
For details see regionImageToEdgeImage in the vigra C++ documentation.
Remove short edges from an edge image.
For details see removeShortEdges in the vigra C++ documentation.
Detect and mark edges in a crack edge image using the Shen/Castan zero-crossing detector.
For details see differenceOfExponentialCrackEdgeImage in the vigra C++ documentation.
Detect and mark edges in an edge image using the Shen/Castan zero-crossing detector.
For details see differenceOfExponentialEdgeImage in the vigra C++ documentation.
watersheds( (object)image [, (int)neighborhood=4 [, (object)seeds=None [, (str)method=’RegionGrowing’ [, (SRGType)terminate=vigra.analysis.SRGType.CompleteGrow [, (int)max_cost=0.0 [, (object)out=None]]]]]]) -> tuple
Compute the watersheds of a 2D image.
- watersheds(image, neighborhood=4, seeds = None, methods = ‘RegionGrowing’,
- terminate=CompleteGrow, threshold=0, out = None) -> (labelimage, max_ragion_label)
Parameters:
- image:
- the image or volume containing the boundary indicator values (high values = high edgeness, dtype=numpy.uint8 or numpy.float32).
- neighborhood:
the pixel neighborhood to be used. Feasible values depend on the dimension and method:
- 2-dimensional data:
- 4 (default) or 8.
- 3-dimensional data:
- 6 (default) or 26
- seeds:
- a label image specifying region seeds, only supported by method ‘RegionGrowing’ (with dtype=numpy.uint32).
- method:
the algorithm to be used for watershed computation. Possible values:
- ‘RegionGrowing’:
- (default) use seededRegionGrowing or seededRegionGrowing3D respectively
- ‘UnionFind:
- use watershedsUnionFind or watersheds3D respectively
- terminate:
when to stop growing. Possible values:
- CompleteGrow:
- (default) grow until all pixels are assigned to a region
- KeepCountours:
- keep a 1-pixel wide contour between all regions, only supported by method ‘RegionGrowing’
- StopAtThreshold:
- stop when the boundary indicator values exceed the threshold given by parameter ‘max_cost’, only supported by method ‘RegionGrowing’
- KeepCountours | StopAtThreshold:
- keep 1-pixel wide contour and stop at given ‘max_cost’, only supported by method ‘RegionGrowing’
- max_cost:
- terminate growing when boundary indicator exceeds this value (ignored when ‘terminate’ is not StopAtThreshold or method is not ‘RegionGrowing’)
- out:
- the label image (with dtype=numpy.uint32) to be filled by the algorithm. It will be allocated by the watershed function if not provided)
The function returns a Python tuple (labelImage, maxRegionLabel)
watersheds( (object)volume [, (int)neighborhood=6 [, (object)seeds=None [, (str)method=’RegionGrowing’ [, (SRGType)terminate=vigra.analysis.SRGType.CompleteGrow [, (int)max_cost=0.0 [, (object)out=None]]]]]]) -> tuple
Compute watersheds of an image using the union find algorithm. If ‘neighborhood’ is ‘None’, it defaults to 8-neighborhood for 2D inputs and 6-neighborhood for 3D inputs.
Calls watersheds() with parameters:
watersheds(image, neighborhood=neighborhood, method='UnionFind', out=out)
The module vigra.geometry contains geometric primitives (such as polygons) and related algorithms.
convexHull( (object)points) -> object
convexHull( (object)points) -> object
Compute the convex hull of a point set.
For details see convexHull in the vigra C++ documentation.
The module vigra.learning will eventually provide a wide range of machine learning tools. Right now, it only contains an implementation of the random forest classifier and probabilistic latent semantic analysis (pLSA) as an example for unsupervised learning.
pLSA( (object)features, (int)nComponents [, (int)nIterations=50 [, (float)minGain=0.0001 [, (bool)normalize=True]]]) -> tuple :
Perform probabilistic latent semantic analysis.
See pLSA in the C++ documentation for detailed information. Note that the feature matrix must have shape (numFeatures * numSamples)!
principleComponents( (object)features, (int)nComponents) -> tuple :
Perform principle component analysis.
See principleComponents in the C++ documentation for detailed information. Note that the feature matrix must have shape (numFeatures * numSamples)!
Constructor:
RandomForest(treeCount = 255, mtry=RF_SQRT, min_split_node_size=1,
training_set_size=0, training_set_proportions=1.0,
sample_with_replacement=True, sample_classes_individually=False,
prepare_online_learning=False)
‘treeCount’ controls the number of trees that are created.
See RandomForest and RandomForestOptions in the C++ documentation for the meaning of the other parameters.
Load from HDF5 file:
RandomForest(filename, pathInFile)
Trains a random Forest using ‘trainData’ and ‘trainLabels’.
and returns the OOB. See the vigra documentation for the meaning af the rest of the paremeters.
Train a random Forest using ‘trainData’ and ‘trainLabels’.
and returns the OOB and the Variable importanceSee the vigra documentation for the meaning af the rest of the paremeters.
Learn online.
Works only if forest has been created with prepare_online_learning=true. Needs the old training data and the new appened, starting at startIndex.
Predict labels on ‘testData’.
The output is an array containing a labels for every test samples.
Predict probabilities for different classes on ‘testData’.
The output is an array containing a probability for every test sample and class.
Re-learn one tree of the forest using ‘trainData’ and ‘trainLabels’.
and returns the OOB. This might be helpful in an online learning setup to improve the classifier.
For more information, refer to RandomForest in the C++ documentation.
Constructor:
RandomForestOld(trainData, trainLabels,
treeCount = 255, mtry=0, min_split_node_size=1,
training_set_size=0, training_set_proportions=1.0,
sample_with_replacement=True, sample_classes_individually=False,)
Construct and train a RandomForest using ‘trainData’ and ‘trainLabels’. ‘treeCount’ controls the number of trees that are created.
See RandomForest and RandomForestOptions in the C++ documentation for the meaning of the other parameters.
The module vigra.noise provides noise estimation and normalization according to a method proposed by Foerstner.
Noise normalization by means of an estimated linear noise model.
For details see linearNoiseNormalization in the vigra C++ documentation.
Determine the noise variance as a function of the image intensity and cluster the results. This operator first calls noiseVarianceEstimation() to obtain a sequence of intensity/variance pairs, which are then clustered using the median cut algorithm. Then the cluster centers (i.e. average variance vs. average intensity) are determined and returned in the result sequence.
Since the length of the resulting array is not known beforhand, it cannot be written into an preallocated array (the “out” argument in most other vigra python functions) . For details see the vigra documentation noiseVarianceClustering.
Determine the noise variance as a function of the image intensity.
Returns an array with the means in the first column and the variances in the second column. Since the length of the resulting array is not known beforhand, it can not be written into an preallocated array (the “out” argument in most other vigra python functions.
For details see the vigra documentation noiseVarianceEstimation.
Noise normalization by means of an estimated non-parametric noise model.
For details see nonparametricNoiseNormalization in the vigra C++ documentation.
Noise normalization by means of an estimated quadratic noise model.
For details see quadraticNoiseNormalization in the vigra C++ documentation.
When you want to write your own vigranumpy extension modules, first make sure that you compile and link with the same versions of numpy and boost_python that your current vigranumpy installation uses. Otherwise, communication between new and existing modules will not work (and even crash). Then follow these steps:
Create the main module source file. This file contains the module’s ‘init’ function. Let’s assume that the module will be called ‘my_module’, and the file is ‘my_module.cxx’. A stub for ‘my_module.cxx’ typically looks like this:
// define PY_ARRAY_UNIQUE_SYMBOL (required by the numpy C-API)
#define PY_ARRAY_UNIQUE_SYMBOL my_module_PyArray_API
// include the vigranumpy C++ API
#include <Python.h>
#include <boost/python.hpp>
#include <vigra/numpy_array.hxx>
#include <vigra/numpy_array_converters.hxx>
... // your includes
... // implementation of your wrapper functions and classes
using namespace boost::python;
// the argument of the init macro must be the module name
BOOST_PYTHON_MODULE_INIT(my_module)
{
// initialize numpy and vigranumpy
vigra::import_vigranumpy();
// export a function
def("my_function", &my_function,
(arg("arg1"), arg("arg2"), ...),
"Documentation");
// export a class and its member functions
class_<MyClass>("MyClass",
"Documentation")
.def("foo", &MyClass::foo,
(arg("arg1"), arg("arg2"), ...),
"Documentation")
;
... // more module functionality (refer to boost_python documentation)
}
When your module uses additional C++ source files, they should start with the following defines:
// this must define the same symbol as the main module file (numpy requirement)
#define PY_ARRAY_UNIQUE_SYMBOL my_module_PyArray_API
#define NO_IMPORT_ARRAY
Implement your wrapper functions. Numpy ndarrays are passed to C++ via the wrapper classes NumpyArray and NumpyAnyArray. You can influence the conversion from Python to C++ by using different instantiations of NumpyArray, as long as the Python array supports the axistags attribute (refer to axis order definitions for the meaning of the term ‘ascending order’):
// We add a 'using' declaration for brevity of our examples.
// In actual code, you should probably prefer explicit namespace qualification.
using namespace vigra;
// Accept any array type and return an arbitrary array type.
// Returning NumpyAnyArray is always safe, because at that point
// C++ no longer cares about the particular type of the array.
NumpyAnyArray foo(NumpyAnyArray array);
// Accept a 3-dimensional float32 array and transpose it
// into ascending axis order ('F' order).
void foo(NumpyArray<3, float> array);
// Accept a 2-dimensional float32 array with an arbitrary number of channels and
// transpose the axes into VIGRA ('V') order (channels are last, other axes ascending).
// Note that the NumpyArray dimension is 3 to account for the channel dimension.
// If the original numpy array has no channel axis, vigranumpy will automatically
// insert a singleton axis.
void foo(NumpyArray<3, Multiband<float> > array);
// Accept a 2-dimensional float32 array that has only a single channel
// (that is, 'array.channels == 1' must hold on the Python side).
// Non-channel axes are transposed into ascending order.
// Note that the NumpyArray dimension is now 2.
void foo(NumpyArray<2, Singleband<float> > array);
// Accept a float32 array that has 2 non-channel dimensions and
// exactly 3 channels (i.e. 'array.channels == 3' on the Python side).
// Non-channel axes are transposed into ascending order.
// Note that the NumpyArray dimension is again 2, but the pixel type is
// now a vector.
// The conversion will only succeed if the channel axis is unstrided on
// the Python side (that is, the following expression is True:
// array.strides[array.channelIndex] == array.dtype.itemsize).
void foo(NumpyArray<2, TinyVector<float, 3> > array);
void foo(NumpyArray<2, RGBValue<float> > array);
Or course, these functions can also be templated.
When your functions return newly allocated arrays, it is usually desirable to transfer the input’s axistags to the output (otherwise, vigranumpy will use defaultAxistags() as a fallback). There is a standard vigranumpy idiom for this task which assumes that the wrapped function has an optional parameter ‘output’ for a possibly pre-allocated output array. The axistags are then transferred by reshaping the output array with a taggedShape() (which is a combination of a shape and axistags):
NumpyAnyArray
foo(NumpyArray<3, Multiband<float32> > input,
NumpyArray<3, Multiband<float32> > output = boost::python::object())
{
// Reshape only if the output array was not explicitly passed in.
// Otherwise, use the output array as is.
output.reshapeIfEmpty(input.taggedShape(),
"error message when shape is unsuitable.");
... // your algorithm
}
It is also possible to modify the tagged shape before it is applied to the output array:
input.taggedShape()
.resize(Shape2(new_width, new_height))
.setChannelCount(new_channel_count)
.setChannelDescription("a description")
The C++ code can be multi-threaded when you unlock Python’s global interpreter lock. After unlocking, your wrapper code must not call any Python functions, so the unlock statement should go after output.reshapeIfEmpty():
NumpyAnyArray
foo(NumpyArray<3, Multiband<float32> > input,
NumpyArray<3, Multiband<float32> > output = boost::python::object())
{
output.reshapeIfEmpty(input.taggedShape(), "Message.");
// Allow parallelization from here on. The destructor of
// _pythread will automatically regain the global interpreter lock
// just before this function returns to Python.
PyAllowThreads _pythread;
... // your algorithm
}
Export your wrapped functions. boost::python::def is called in its usual way, with one simple extension: Since vigranumpy does not know which NumpyArray variants you are going to use, appropriate converter functions between Python and C++ must be registered on demand. You do this by enclosing your function pointer into a call to the ‘registerConverters()’ function:
// in the module's init function
def("my_function", vigra::registerConverters(&my_function),
(arg("arg1"), ...),
"Documentation");
If you need more information, it is always a good idea to look at the source code of the existing vigranumpy modules.