[ VIGRA Homepage | Function Index | Class Index | Namespaces | File List | Main Page ]

details GridGraph< N, DirectedTag > Class Template Reference VIGRA

Define a grid graph in arbitrary dimensions. More...

#include <vigra/multi_gridgraph.hxx>

Inherits GridGraphBase< N, DirectedTag >.

Classes

class  ArcMap
 Type of an arc property map that maps arc descriptor objects onto property values of type T (API: LEMON). More...
 
class  EdgeMap
 Type of an edge property map that maps edge descriptor objects onto property values of type T (API: LEMON). More...
 
class  InDegMap
 Type of a property map that returns the number of incoming edges of a given node (API: LEMON, use via lemon::InDegMap<Graph>). More...
 
class  IndexMap
 Type of a property map that returns the coordinate of a given node (API: LEMON). More...
 
class  NodeMap
 Type of a node property map that maps node descriptor objects onto property values of type T (API: LEMON). More...
 
class  OutDegMap
 Type of a property map that returns the number of outgoing edges of a given node (API: LEMON, use via lemon::OutDegMap<Graph>). More...
 
struct  traversal_category
 Define several graph tags related to graph traversal (API: boost::graph, use via boost::graph_traits<Graph>::traversal_category). More...
 

Public Types

typedef
GridGraphNeighborIterator< N > 
adjacency_iterator
 Iterator over the neighbor vertices of a given vertex (API: boost::graph, use via boost::graph_traits<Graph>::adjacency_iterator).
 
typedef GridGraphArcDescriptor< N > Arc
 The arc (directed edge) descriptor (API: LEMON).
 
typedef GridGraphArcIterator
< N, false > 
ArcIt
 Iterator over the acrs (directed edges) of a node (API: LEMON).
 
typedef
GridGraphNeighborIterator< N,
true > 
back_neighbor_vertex_iterator
 Iterator over only those neighbor vertices of a given vertex that have smaller ID (API: VIGRA).
 
typedef MultiArrayIndex degree_size_type
 Type to specify number of neighbors (API: boost::graph, use via boost::graph_traits<Graph>::degree_size_type).
 
typedef DirectedTag directed_category
 Is the graph directed or undirected ? (API: boost::graph, use via boost::graph_traits<Graph>::directed_category).
 
typedef MultiArrayShape< N+1 >
::type 
Edge
 The edge descriptor (API: LEMON).
 
typedef GridGraphArcDescriptor< N > edge_descriptor
 The edge descriptor (API: boost::graph, use via boost::graph_traits<Graph>::edge_descriptor).
 
typedef GridGraphArcIterator
< N,!is_directed
edge_iterator
 Iterator over the edges of a graph (API: boost::graph, use via boost::graph_traits<Graph>::edge_iterator).
 
typedef
boost_graph::disallow_parallel_edge_tag 
edge_parallel_category
 The graph does not allow multiple edges between the same vertices (API: boost::graph, use via boost::graph_traits<Graph>::edge_parallel_category).
 
typedef MultiArrayShape< N+1 >
::type 
edge_propmap_shape_type
 Shape type of an edge property map (must have one additional dimension).
 
typedef GridGraphEdgeIterator
< N,!is_directed
EdgeIt
 Iterator over the edges of the graph (API: LEMON).
 
typedef MultiArrayIndex edges_size_type
 Type to specify number of edges (API: boost::graph, use via boost::graph_traits<Graph>::edges_size_type).
 
typedef GridGraphInArcIterator< N > in_edge_iterator
 Iterator over the incoming edges of a given vertex (API: boost::graph, use via boost::graph_traits<Graph>::in_edge_iterator).
 
typedef GridGraphInArcIterator< N > InArcIt
 Iterator over the incoming arcs of a node (API: LEMON).
 
typedef
GridGraphOutEdgeIterator< N,
true > 
IncBackEdgeIt
 Iterator over only those incident edges of a node that end in a node with smaller ID (API: VIGRA).
 
typedef
GridGraphOutEdgeIterator< N > 
IncEdgeIt
 Iterator over the incident edges of a node (API: LEMON).
 
typedef MultiArrayIndex index_type
 Type of node and edge IDs.
 
typedef adjacency_iterator neighbor_vertex_iterator
 Same as adjacency_iterator (API: VIGRA).
 
typedef vertex_descriptor Node
 The Node descriptor (API: LEMON).
 
typedef vertex_iterator NodeIt
 Iterator over all nodes of the graph (API: LEMON).
 
typedef
GridGraphOutArcIterator< N,
true > 
out_back_edge_iterator
 Iterator over only those outgoing edges of a given vertex that go to vertices with smaller IDs (API: VIGRA).
 
typedef
GridGraphOutArcIterator< N > 
out_edge_iterator
 Iterator over the outgoing edges of a given vertex (API: boost::graph, use via boost::graph_traits<Graph>::out_edge_iterator).
 
typedef
GridGraphOutArcIterator< N > 
OutArcIt
 Iterator over the outgoing edges of a node (API: LEMON).
 
typedef
GridGraphOutArcIterator< N,
true > 
OutBackArcIt
 Iterator over only those outgoing edges of a node that end in a node with smaller ID (API: VIGRA).
 
typedef MultiArrayShape< N >::type shape_type
 Shape type of the graph and a node property map.
 
typedef shape_type vertex_descriptor
 The vertex descriptor (API: boost::graph, use via boost::graph_traits<Graph>::vertex_descriptor).
 
typedef
MultiCoordinateIterator< N > 
vertex_iterator
 Iterator over the vertices of the graph (API: boost::graph, use via boost::graph_traits<Graph>::vertex_iterator).
 
typedef boost_graph::no_property vertex_property_type
 The graph does not define internal property maps (API: boost::graph, use via boost::graph_traits<Graph>::vertex_property_type).
 
typedef MultiArrayIndex vertices_size_type
 Type to specify number of vertices (API: boost::graph, use via boost::graph_traits<Graph>::vertices_size_type).
 

Public Member Functions

Arc arcFromId (index_type i) const
 Get an arc descriptor for the given arc ID i (API: LEMON). More...
 
edges_size_type arcNum () const
 Get the number of arc in this graph (API: LEMON).
 
degree_size_type back_degree (vertex_descriptor const &v) const
 Get the number of outgoing backward edges of the given vertex (API: VIGRA).
 
Node baseNode (IncEdgeIt const &e) const
 Return the start node of the edge the given iterator is referring to (API: LEMON).
 
Node baseNode (IncBackEdgeIt const &e) const
 Return the start node of the edge the given iterator is referring to (API: VIGRA).
 
Node baseNode (OutArcIt const &a) const
 Return the start node of the edge the given iterator is referring to (API: LEMON).
 
Node baseNode (OutBackArcIt const &a) const
 Return the start node of the edge the given iterator is referring to (API: VIGRA).
 
degree_size_type degree (vertex_descriptor const &v) const
 Get the total number of edges (incoming plus outgoing) of the given vertex (convenience function,
the boost::graph API provides the free function boost::degree(v, graph),
LEMON has no analogue).
 
Arc direct (Edge const &e, bool forward) const
 Create an arc for the given edge e, oriented along the edge's natural (forward = true) or reversed (forward = false) direction (API: LEMON).
 
Arc direct (Edge const &e, Node const &n) const
 Create an arc for the given edge e oriented so that node n is the starting node of the arc (API: LEMON), or return lemon::INVALID if the edge is not incident to this node.
 
bool direction (Arc const &a) const
 Return true when the arc is looking on the underlying edge in its natural (i.e. forward) direction, false otherwise (API: LEMON).
 
std::pair< edge_descriptor, bool > edge (vertex_descriptor const &u, vertex_descriptor const &v) const
 Get a descriptor for the edge connecting vertices u and v,
or (lemon::INVALID, false) if no such edge exists (convenience function,
the boost::graph API provides the free function boost::edge(u, v, graph)).
 
Edge edgeFromId (index_type i) const
 Get the edge descriptor for the given edge ID i (API: LEMON). More...
 
edges_size_type edgeNum () const
 Get the number of edges in this graph (API: LEMON).
 
Arc findArc (Node const &u, Node const &v, Arc const &=lemon::INVALID) const
 Get a descriptor for the arc connecting vertices u and v,
or lemon::INVALID if no such edge exists (API: LEMON).
 
Edge findEdge (Node const &u, Node const &v, Edge const &=lemon::INVALID) const
 Get a descriptor for the edge connecting vertices u and v,
or lemon::INVALID if no such edge exists (API: LEMON).
 
degree_size_type forward_degree (vertex_descriptor const &v) const
 Get the number of outgoing forward edges of the given vertex (API: VIGRA).
 
back_neighbor_vertex_iterator get_back_neighbor_vertex_end_iterator (vertex_descriptor const &v) const
 Get an iterator pointing beyond the range of backward neighbors of the given vertex (API: VIGRA,
in analogy to the boost::graph API, we also provide a free function boost::back_adjacent_vertices(v, g),
and LEMON just uses lemon::INVALID instead).
 
back_neighbor_vertex_iterator get_back_neighbor_vertex_iterator (vertex_descriptor const &v) const
 Get an iterator pointing to the first backward neighbor of the given vertex (API: VIGRA,
in analogy to the boost::graph API, we also provide a free function boost::back_adjacent_vertices(v, g),
and the LEMON analogue is Graph::OutBackArcIt(graph, v)).
 
edge_iterator get_edge_end_iterator () const
 Get edge iterator pointing beyond the valid range of edges of this graph (convenience function,
the boost::graph API provides the free function boost::vertices(graph),
LEMON uses the special value lemon::INVALID instead).
 
edge_iterator get_edge_iterator () const
 Get edge iterator pointing to the first edge of the graph (convenience function,
the boost::graph API provides the free function boost::edges(graph),
LEMON uses Graph::EdgeIt(graph)).
 
in_edge_iterator get_in_edge_end_iterator (vertex_descriptor const &v) const
 Get an iterator pointing beyond the range of incoming edges of the given vertex (convenience function,
the boost::graph API provides the free function boost::in_edges(v, graph),
LEMON uses the special value lemon::INVALID instead).
 
in_edge_iterator get_in_edge_iterator (vertex_descriptor const &v) const
 Get an iterator pointing to the first incoming edge of the given vertex (convenience function,
the boost::graph API provides the free function boost::in_edges(v, graph),
LEMON uses Graph::InArcIt(g, v)).
 
neighbor_vertex_iterator get_neighbor_vertex_end_iterator (vertex_descriptor const &v) const
 Get an iterator pointing beyond the range of neighbors of the given vertex (convenience function,
the boost::graph API provides the free function boost::adjacent_vertices(v, graph),
LEMON uses the speical value lemon::INVALID instead).
 
neighbor_vertex_iterator get_neighbor_vertex_iterator (vertex_descriptor const &v) const
 Get an iterator pointing to the first neighbor of the given vertex (convenience function,
the boost::graph API provides the free function boost::adjacent_vertices(v, graph),
LEMON uses Graph::ArcIt(g, v)).
 
out_back_edge_iterator get_out_back_edge_end_iterator (vertex_descriptor const &v) const
 Get an iterator pointing beyond the range of outgoing backward edges of the given vertex (API: VIGRA,
in analogy to the boost::graph API, we also provide a free function boost::out_back_edges(v, g),
and LEMON uses the special value lemon::INVALID instead).
 
out_back_edge_iterator get_out_back_edge_iterator (vertex_descriptor const &v) const
 Get an iterator pointing to the first outgoing backward edge of the given vertex (API: VIGRA,
in analogy to the boost::graph API, we also provide a free function boost::out_back_edges(v, g),
and the LEMON analogue is Graph::IncBackEdgeIt(graph, v)).
 
out_edge_iterator get_out_edge_end_iterator (vertex_descriptor const &v) const
 Get an iterator pointing beyond the range of outgoing edges of the given vertex (convenience function,
the boost::graph API provides the free function boost::out_edges(v, graph),
LEMON uses the special value lemon::INVALID instead).
 
out_edge_iterator get_out_edge_iterator (vertex_descriptor const &v) const
 Get an iterator pointing to the first outgoing edge of the given vertex (convenience function,
the boost::graph API provides the free function boost::out_edges(v, graph),
LEMON uses Graph::OutArcIt(g, v)).
 
vertex_iterator get_vertex_end_iterator () const
 Get vertex iterator pointing beyond the valid range of vertices of this graph (convenience function,
the boost::graph API provides the free function boost::vertices(graph),
LEMON uses the special value lemon::INVALID instead).
 
vertex_iterator get_vertex_iterator () const
 Get vertex iterator pointing to the first vertex of this graph (convenience function,
the boost::graph API provides the free function boost::vertices(graph),
LEMON uses Graph::NodeIt(graph)).
 
vertex_iterator get_vertex_iterator (vertex_descriptor const &v) const
 Get vertex iterator pointing to the given vertex (API: VIGRA).
 
 GridGraph (shape_type const &shape, NeighborhoodType ntype=DirectNeighborhood)
 Construct a grid graph with given shape and neighborhood type ntype. More...
 
index_type id (Node const &v) const
 Get the ID (i.e. scan-order index) for node desciptor v (API: LEMON).
 
index_type id (Edge const &e) const
 Get the ID (i.e. scan-order index in an edge property map) for the given edges descriptor e (API: LEMON).
 
index_type id (Arc const &a) const
 Get the ID (i.e. scan-order index an an arc property map) for the given ar a (API: LEMON).
 
degree_size_type in_degree (vertex_descriptor const &v) const
 Get the number of incoming edges of the given vertex (convenience function,
the boost::graph API provides the free function boost::in_degree(v, graph),
LEMON uses a special property map lemon::InDegMap<Graph>).
 
IndexMap indexMap () const
 Create a property map that returns the coordinate of each node (API: LEMON GridGraph).
 
index_type maxArcId () const
 Get the maximal ID af any arc in this graph (API: LEMON).
 
index_type maxEdgeId () const
 Get the maximum ID of any edge in this graph (API: LEMON).
 
index_type maxNodeId () const
 Get the maximum ID of any node in this graph (API: LEMON).
 
Node nodeFromId (index_type i) const
 Get node descriptor for given node ID i (API: LEMON). More...
 
vertices_size_type nodeNum () const
 Get the number of nodes in this graph (API: LEMON).
 
edges_size_type num_edges () const
 Get the number of edges in this graph (convenience function, boost::graph API provides the free function boost::num_edges(graph)).
 
vertices_size_type num_vertices () const
 Get the number of vertices in this graph (convenience function, the boost::graph API provides the free function boost::num_vertices(graph)).
 
Arc oppositeArc (Arc const &a) const
 Create an arc referring to the same edge as the given arc a, but with reversed direction (API: LEMON).
 
Node oppositeNode (Node const &n, Edge const &e) const
 Return the opposite node of the given node n along edge e (API: LEMON), or return lemon::INVALID if the edge is not incident to this node.
 
degree_size_type out_degree (vertex_descriptor const &v) const
 Get the number of outgoing edges of the given vertex (convenience function,
the boost::graph API provides the free function boost::out_degree(v, graph),
LEMON uses a special property map lemon::OutDegMap<Graph>).
 
Node const & pos (Node const &v) const
 Get the grid cordinate of the given node v (convenience function).
 
Node runningNode (IncEdgeIt const &e) const
 Return the end node of the edge the given iterator is referring to (API: LEMON).
 
Node runningNode (IncBackEdgeIt const &e) const
 Return the end node of the edge the given iterator is referring to (API: VIGRA).
 
Node runningNode (OutArcIt const &a) const
 Return the end node of the edge the given iterator is referring to (API: LEMON).
 
Node runningNode (OutBackArcIt const &a) const
 Return the end node of the edge the given iterator is referring to (API: VIGRA).
 
Node source (Arc const &a) const
 Get the start node of the given arc a (API: LEMON).
 
Node target (Arc const &a) const
 Get the end node of the given arc a (API: LEMON).
 
Node u (Edge const &e) const
 Get the start node of the given edge e (API: LEMON,
the boost::graph API provides the free function boost::source(e, graph)).
 
Node v (Edge const &e) const
 Get the end node of the given edge e (API: LEMON,
the boost::graph API provides the free function boost::target(e, graph)).
 

Static Public Attributes

static const unsigned int dimension = N
 Dimension of the grid.
 
static const bool is_directed = IsSameType<DirectedTag, directed_tag>::value
 'true' if the graph is directed (API: boost::graph)
 

Detailed Description

template<unsigned int N, class DirectedTag>
class vigra::GridGraph< N, DirectedTag >

Define a grid graph in arbitrary dimensions.

A GridGraph connects each pixel to its direct or indirect neighbors. Direct neighbors are the adjacent point along the coordinate axes, whereas indirect neighbors include the diagonally adjacent points. Thus, direct neighbors correspond to the 4-neighborhood in 2D and the 6-neighborhood in 3D, whereas indirect neighbors correspond to the 8-neighborhood and 26-neighborhood respectively. The GridGraph class extends this scheme to arbitrary dimensions. While the dimension must be defined at compile time via the template parameter N, the desired neighborhood can be chosen at runtime in the GridGraph's constructor. The shape of the grid is also specified at runtime in terms of a suitable shape object.

Another choice to be made at compile time is whether the graph should be directed or undirected. This is defined via the DirectedTag template parameter which can take the values directed_tag or undirected_tag (default).

The main difficulty in a grid graph is to skip the missing neighbors of the pixels near the grid's border. For example, the upper left pixel in a 2D grid has only two direct neighbors instead of the usual four. The GridGraph class uses a precomputed set of internal look-up tables to efficiently determine the appropriate number and location of the existing neighbors. A key design decision to make this fast was to give up the requirement that edge IDs are contiguous integers (as in LEMON's implementation of a 2D grid graph, which became very complicated and slow by enforcing this requirement). Instead, edges are numbered as if all nodes (including nodes at the grid's border) had the same degree. Since edges crossing the border are not actually present in the graph, these IDs are missing, leading to gaps in the sequence of IDs.

For the sake of compatibility, the GridGraph class implements two popular graph APIs: the one defined by the boost::graph library and the one defined by the LEMON library. Their basic principles are very similar: The graph's nodes, edges and adjacency structure are accessed via a set of iterators, whereas additional properties (like node and edge weights) are provided via property maps that are indexed with suitable node or edge descriptor objects returned by the iterators.

Specifically, GridGraph implements the requirements of the following concepts of the boost::graph API: Graph, IncidenceGraph, BidirectionalGraph, AdjacencyGraph, VertexListGraph, EdgeListGraph, and AdjacencyMatrix. Likewise, it supports the concepts Graph and Digraph of the LEMON API. The property maps associated with a GridGraph support the boost concepts ReadablePropertyMap, WritablePropertyMap, ReadWritePropertyMap, and LvaluePropertyMap as well as the LEMON concepts ReadMap, WriteMap, ReadWriteMap, and ReferenceMap.

VIGRA's GridGraph class is designed such that multi-dimensional coordinates (i.e. TinyVector<MultiArrayIndex>) serve as descriptor objects, which means that normal MultiArrays or MultiArrayViews can serve as property maps in most situations. Thus, node properties like a foreground probability for foreground/background segmentation can simply be stored as normal images.

Since the boost::graph and LEMON APIs differ in how they call corresponding functionality (e.g., they use the terms vertex and node respectively in an exactly synonymous way), most GridGraph helper classes and functions are exposed under two different names. To implement your own algorithms, you can choose the API you like most. VIGRA adopts the convention that algorithms using the boost::graph API go into the namespace vigra::boost_graph, while those using the LEMON API are placed into the namespace vigra::lemon_graph. This helps to avoid name clashes when the two APIs happen to use the same name for different things. The documentation of the GridGraph members specifies which API the respective function or class belongs to. Please consult the documentation of these libraries for tutorial introductions and full reference of the respective APIs.

VIGRA adds an important new concept of its own: the back neighborhood and associated adjacency and edge iterators. The back neighborhood of a given vertex with ID i is the set of all neighbor vertices whose ID is smaller than i. This concept is useful if you work on the grid graph's vertices in scan order and want to access just those neighbors that have already been processed. Connected components labeling is a typical example of an algorithm that can take advantage of this possibility. In principle, a back neighborhood iterator could be implemented as a filter iterator that simply skips all neighbors with inappropriate IDs. However, in case of grid graphs it is more efficient to provide a special iterator for this purpose.

Usage:

#include <vigra/multi_gridgraph.hxx>
Namespace: vigra

At present, the GridGraph class is mainly used internally to implement image analysis functions for arbitrary dimensional arrays (e.g. detection of local extrema, connected components labeling, watersheds, SLIC superpixels). For example, a dimension-independent algorithm to detect local maxima using the LEMON API might look like this:

namespace vigra { namespace lemon_graph {
template <class Graph, class InputMap, class OutputMap>
void
localMaxima(Graph const & g,
InputMap const & src,
OutputMap & local_maxima,
typename OutputMap::value_type marker)
{
// define abreviations for the required iterators
typedef typename Graph::NodeIt graph_scanner;
typedef typename Graph::OutArcIt neighbor_iterator;
// iterate over all nodes (i.e. pixels)
for (graph_scanner node(g); node != INVALID; ++node)
{
// remember the value of the current node
typename InputMap::value_type current = src[*node];
// iterate over all neighbors of the current node
bool is_local_maximum = true;
for (neighbor_iterator arc(g, node); arc != INVALID; ++arc)
{
// if a neighbor has larger value, the center node is not a local maximum
if (current < src[g.target(*arc)])
is_local_maximum = false;
}
// mark the node when it is a local maximum
if (is_local_maximum)
local_maxima[*node] = marker;
}
}
}} // namespace vigra::lemon_graph

It should be noted that this algorithm will work for any LEMON-compatible graph class, not just our GridGraph. When implemented in terms of the boost::graph API, the same algorithm looks like this:

namespace vigra { namespace boost_graph {
template <class Graph, class InputMap, class OutputMap>
void
localMaxima(Graph const & g,
InputMap const & src,
OutputMap & local_maxima,
typename property_traits<OutputMap>::value_type marker)
{
// define abreviations and variables for the required iterators
typedef typename graph_traits<Graph>::vertex_iterator graph_scanner;
typedef typename graph_traits<Graph>::adjacency_iterator neighbor_iterator;
graph_scanner node, node_end;
neighbor_iterator arc, arc_end;
// iterate over all nodes (i.e. pixels)
tie(node, node_end) = vertices(g);
for (; node != node_end; ++node)
{
// remember the value of the current node
typename property_traits<InputMap>::value_type current = get(src, *node);
// iterate over all neighbors of the current node
bool is_local_maximum = true;
tie(arc, arc_end) = adjacent_vertices(*node, g);
for (;arc != arc_end; ++arc)
{
// if a neighbor has larger value, the center node is not a local maximum
if (current < get(src, *arc))
is_local_maximum = false;
}
// mark the node when it is a local maximum
if (is_local_maximum)
put(local_maxima, *node, marker);
}
}
}} // namespace vigra::boost_graph

It can be seen that the differences between the APIs are mainly syntactic (especially note that boost::graph users traits classes and free functions, whereas LEMON uses nested typedefs and member functions). Either of these algorithms can now serve as the backend of a local maxima detector for MultiArrayViews:

namespace vigra {
template <unsigned int N, class T1,
class T2>
void
localMaxima(MultiArrayView<N, T1> const & src,
MultiArrayView<N, T2> local_maxima,
T2 marker,
{
vigra_precondition(src.shape() == local_maxima.shape(),
"localMinMax(): shape mismatch between input and output.");
// create a grid graph with appropriate shape and desired neighborhood
GridGraph<N, undirected_tag> graph(src.shape(), neighborhood);
// forward the call to the graph-based algorithm, using
// the given MultiArrayViews as property maps
lemon_graph::localMaxima(graph, src, local_maxima, marker);
}
} // namespace vigra

A slightly enhanced version of this code is actually used to implement this functionality in VIGRA.

Examples:
graph_agglomerative_clustering.cxx.

Constructor & Destructor Documentation

GridGraph ( shape_type const &  shape,
NeighborhoodType  ntype = DirectNeighborhood 
)

Construct a grid graph with given shape and neighborhood type ntype.

The shape must have type MultiArrayShape<N>::type with the appropriate dimension N. The neighborhood type can take the values DirectNeighborhood to use only the axis-aligned edges (2N-neighborhood) and IndirectNeighborhood to use all diagonal edges as well ((3N-1)-neighborhood).

Member Function Documentation

Node nodeFromId ( index_type  i) const

Get node descriptor for given node ID i (API: LEMON).

Return Node(lemon::INVALID) when the ID does not exist in this graph.

Edge edgeFromId ( index_type  i) const

Get the edge descriptor for the given edge ID i (API: LEMON).

Return Edge(lemon::INVALID) when the ID does not exist in this graph.

Arc arcFromId ( index_type  i) const

Get an arc descriptor for the given arc ID i (API: LEMON).

Return Arc(lemon::INVALID) when the ID does not exist in this graph.


The documentation for this class was generated from the following files:

© Ullrich Köthe (ullrich.koethe@iwr.uni-heidelberg.de)
Heidelberg Collaboratory for Image Processing, University of Heidelberg, Germany

html generated using doxygen and Python
vigra 1.11.1 (Fri May 19 2017)