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

details vigra Namespace Reference VIGRA

Namespaces

 acc
 Efficient computation of object statistics.
 
 linalg
 Linear algebra functions.
 
 multi_math
 
 metrics
 Define functors for various common metrics.
 
 FourNeighborhood
 
 EightNeighborhood
 
 rf3
 Random forest version 3.
 
 Neighborhood3DSix
 
 Neighborhood3DTwentySix
 

Classes

class  AdjacencyListGraph
 undirected adjacency list graph in the LEMON API More...
 
class  AffineMotionEstimationOptions
 Option object for affine registration functions. More...
 
class  Any
 Typesafe storage of arbitrary values. More...
 
class  ArrayOfRegionStatistics
 Calculate statistics for all regions of a labeled image. More...
 
class  ArrayVector
 
class  ArrayVectorView
 
class  BasicImage
 Fundamental class template for images. More...
 
class  BasicImageIterator
 
class  BasicImageIteratorBase
 
class  BasicImageView
 BasicImage using foreign memory. More...
 
class  BestGiniOfColumn
 
class  BilinearInterpolatingAccessor
 Bilinear interpolation at non-integer positions. More...
 
class  BinaryForest
 BinaryForest stores a collection of rooted binary trees. More...
 
class  BlockwiseConvolutionOptions
 
class  BlockwiseLabelOptions
 
class  BlockwiseOptions
 
class  BlueAccessor
 
class  Box
 Represent an n-dimensional box as a (begin, end) pair. Depending on the value type, end() is considered to be outside the box (as in the STL, for integer types), or inside (for floating point types). size() will always be end() - begin(). More...
 
class  BrightnessContrastFunctor
 Adjust brightness and contrast of an image. More...
 
class  BSpline
 
class  BSplineBase
 
class  BucketQueue
 Priority queue implemented using bucket sort. More...
 
class  CatmullRomSpline
 
class  ChangeablePriorityQueue
 Heap-based changable priority queue with a maximum number of elemements. More...
 
class  ChunkedArray
 Interface and base class for chunked arrays. More...
 
class  ChunkedArrayCompressed
 
class  ChunkedArrayFull
 
class  ChunkedArrayHDF5
 
class  ChunkedArrayLazy
 
class  ChunkedArrayOptions
 Option object for ChunkedArray construction. More...
 
struct  ChunkedArrayTag
 
class  ChunkedArrayTmpFile
 
class  ClusteringOptions
 Options object for hierarchical clustering. More...
 
class  ColumnIterator
 Iterator adapter to linearly access columns. More...
 
class  ConstBasicImageIterator
 
class  ConstImageIterator
 Standard 2D random access const iterator for images that store the data as a linear array. More...
 
class  ConstStridedImageIterator
 Const iterator to be used when pixels are to be skipped. More...
 
class  ConstValueIterator
 Iterator that always returns the constant specified in the constructor. More...
 
class  ConvexPolytope
 
class  ConvolutionOptions
 Options class template for convolutions. More...
 
class  CoordinateConstValueAccessor
 Forward accessor to the value() part of the values an iterator points to. More...
 
class  CoscotFunction
 
class  CountingIterator
 Iterator that counts upwards or downwards with a given step size. More...
 
struct  CoupledArrays
 
class  CoupledHandle
 
struct  CoupledIteratorType
 
class  CoupledScanOrderIterator
 Iterate over multiple images simultaneously in scan order. More...
 
class  CrackContourCirculator
 Circulator that walks around a given region. More...
 
class  Diff2D
 Two dimensional difference vector. More...
 
class  DiffusivityFunctor
 Diffusivity functor for non-linear diffusion. More...
 
class  Dist2D
 
struct  DistancePowerFunctor
 
class  DT_StackEntry
 
class  EarlyStoppStd
 Standard early stopping criterion. More...
 
class  Edgel
 
class  EnhancedFrostFunctor
 This function tries to reduce the speckle noise of an image by applying the Enhanced Frost filter. More...
 
class  EnhancedLeeFunctor
 This function tries to reduce the speckle noise of an image by applying the Enhanced Lee filter. More...
 
class  EntropyCriterion
 
class  FFTWComplex
 Wrapper class for the FFTW complex types 'fftw_complex'. More...
 
class  FFTWConvolvePlan
 
class  FFTWCorrelatePlan
 
class  FFTWImaginaryAccessor
 
class  FFTWLogMagnitudeAccessor
 
class  FFTWMagnitudeAccessor
 
class  FFTWPhaseAccessor
 
class  FFTWPlan
 
class  FFTWRealAccessor
 
class  FFTWSquaredMagnitudeAccessor
 
class  FFTWWriteRealAccessor
 
class  FilterIterator
 This iterator creates a view of another iterator and skips elements that do not fulfill a given predicate. More...
 
class  FindAverage
 Find the average pixel value in an image or ROI. More...
 
class  FindAverageAndVariance
 Find the average pixel value and its variance in an image or ROI. More...
 
class  FindBoundingRectangle
 Calculate the bounding rectangle of an ROI in an image. More...
 
class  FindMinMax
 Find the minimum and maximum pixel value in an image or ROI. More...
 
class  FindROISize
 Calculate the size of an ROI in an image. More...
 
class  FindSum
 Find the sum of the pixel values in an image or ROI. More...
 
class  FixedPoint
 
class  FixedPoint16
 
class  FunctorTraits
 Export associated information for a functor. More...
 
class  GaborFilterFamily
 Family of gabor filters of different scale and direction. More...
 
class  GammaFunctor
 Perform gamma correction of an image. More...
 
class  GammaMAPFunctor
 This function tries to reduce the speckle noise of an image by applying the Gamma Maximum A Posteriori (MAP) filter. More...
 
class  Gaussian
 
class  GiniCriterion
 
class  GrayToRGBAccessor
 
class  GreenAccessor
 
class  GridGraph
 Define a grid graph in arbitrary dimensions. More...
 
class  HDF5DisableErrorOutput
 Temporarily disable HDF5's native error output. More...
 
class  HDF5File
 Access to HDF5 files. More...
 
class  HDF5Handle
 Wrapper for unique hid_t objects. More...
 
class  HDF5HandleShared
 Wrapper for shared hid_t objects. More...
 
class  HDF5ImportInfo
 Argument object for the function readHDF5(). More...
 
struct  HierarchicalIteratorType
 
class  HistogramOptions
 Set histogram options. More...
 
class  ImageArray
 Fundamental class template for arrays of equal-sized images. More...
 
class  ImageExportInfo
 Argument object for the function exportImage(). More...
 
class  ImageImportInfo
 Argument object for the function importImage(). More...
 
class  ImageIterator
 Standard 2D random access iterator for images that store the data in a linear array. More...
 
class  ImageIteratorBase
 Base class for 2D random access iterators. More...
 
class  ImagePyramid
 Class template for logarithmically tapering image pyramids. More...
 
class  IteratorAdaptor
 Quickly create 1-dimensional iterator adapters. More...
 
struct  IteratorTraits
 Export associated information for each image iterator. More...
 
class  Kernel1D
 Generic 1 dimensional convolution kernel. More...
 
class  Kernel2D
 Generic 2 dimensional convolution kernel. More...
 
class  KuanFunctor
 This function tries to reduce the speckle noise of an image by applying the Kuan filter. More...
 
class  Lab2RGBFunctor
 Convert perceptual uniform CIE L*a*b* into linear (raw) RGB. More...
 
class  Lab2RGBPrimeFunctor
 Convert perceptual uniform CIE L*a*b* into non-linear (gamma corrected) R'G'B'. More...
 
class  Lab2XYZFunctor
 Convert perceptual uniform CIE L*a*b* into standardized tri-stimulus XYZ. More...
 
class  LabelOptions
 Option object for labelMultiArray(). More...
 
class  LastValueFunctor
 Stores and returns the last value it has seen. More...
 
class  LeeFunctor
 This function tries to reduce the speckle noise of an image by applying the basic Lee filter. More...
 
class  LineIterator
 Iterator adapter to iterate along an arbitrary line on the image. More...
 
class  LocalMinmaxOptions
 Options object for localMinima() and localMaxima(). More...
 
class  Luv2RGBFunctor
 Convert perceptual uniform CIE L*u*v* into linear (raw) RGB. More...
 
class  Luv2RGBPrimeFunctor
 Convert perceptual uniform CIE L*u*v* into non-linear (gamma corrected) R'G'B'. More...
 
class  Luv2XYZFunctor
 Convert perceptual uniform CIE L*u*v* into standardized tri-stimulus XYZ. More...
 
class  MagnitudeFunctor
 
class  MappedBucketQueue
 Priority queue implemented using bucket sort (STL compatible). More...
 
class  MergeGraphAdaptor
 undirected graph adaptor for edge contraction and feature merging More...
 
struct  MeshGridAccessor
 
class  MultiArray
 Main MultiArray class containing the memory management. More...
 
class  MultiArrayNavigator
 A navigator that provides access to the 1D subranges of an n-dimensional range given by a vigra::MultiIterator and an nD shape. More...
 
class  MultiArrayShape
 
class  MultiArrayView
 Base class for, and view to, vigra::MultiArray. More...
 
class  MultiBlocking
 
class  MultiCoordinateIterator
 Iterate over a virtual array where each element contains its coordinate. More...
 
class  MultiCoordinateNavigator
 A navigator that provides access to the 1D subranges of an n-dimensional range given by an nD shape. More...
 
class  MultiImageAccessor2
 Access two images simultaneously. More...
 
class  MultiIterator
 A multi-dimensional hierarchical iterator to be used with vigra::MultiArrayView if it is not strided. More...
 
class  NeighborhoodCirculator
 Circulator that walks around a given location in a given image. More...
 
class  NeighborOffsetCirculator
 Circulator that walks around a given location. More...
 
class  Node< e_ConstProbNode >
 
class  NodeBase
 
class  NoiseNormalizationOptions
 Pass options to one of the noise normalization functions. More...
 
class  NonlinearLSQOptions
 Pass options to nonlinearLeastSquares(). More...
 
class  NormalRandomFunctor
 
class  NumpyAnyArray
 
class  NumpyArray
 
class  ParallelOptions
 Option base class for parallel algorithms. More...
 
class  PLSAOptions
 Option object for the pLSA algorithm. More...
 
class  Point2D
 Two dimensional point or position. More...
 
class  Polygon
 
class  Polynomial
 
class  PolynomialView
 
class  Polytope
 Represent an n-dimensional polytope. More...
 
class  PriorityQueue
 Heap-based priority queue compatible to BucketQueue. More...
 
class  ProblemSpec
 problem specification class for the random forest. More...
 
class  Processor
 
class  Processor< ClassificationTag, LabelType, T1, C1, T2, C2 >
 
class  Processor< RegressionTag, LabelType, T1, C1, T2, C2 >
 
class  PropertyMap
 The PropertyMap is used to store Node or Arc information of graphs. More...
 
class  PropertyMap< KEYTYPE, MAPPEDTYPE, IndexVectorTag >
 Specialization of PropertyMap that stores the elements in a vector (size = number of stored elements). An additional index vector is needed for bookkeeping (size = max node id of stored elements). More...
 
class  PropertyMap< KEYTYPE, MAPPEDTYPE, VectorTag >
 Specialization of PropertyMap that stores the elements in a vector (size = max node id of stored elements). More...
 
class  Quaternion
 
class  RandomForest
 Random forest version 2 (see also vigra::rf3::RandomForest for version 3) More...
 
class  RandomForestClassCounter
 
class  RandomForestOptions
 Options object for the random forest. More...
 
class  RandomNumberGenerator
 
class  Rational
 
class  Rect2D
 Two dimensional rectangle. More...
 
class  RedAccessor
 
class  ReduceFunctor
 Apply a functor to reduce the dimensionality of an array. More...
 
class  RestrictedNeighborhoodCirculator
 Circulator that walks around a given location in a given image, using a restricted neighborhood. More...
 
class  RGB2LabFunctor
 Convert linear (raw) RGB into perceptual uniform CIE L*a*b*. More...
 
class  RGB2LuvFunctor
 Convert linear (raw) RGB into perceptual uniform CIE L*u*v*. More...
 
class  RGB2RGBPrimeFunctor
 Convert linear (raw) RGB into non-linear (gamma corrected) R'G'B'. More...
 
class  RGB2sRGBFunctor
 Convert linear (raw) RGB into standardized sRGB. More...
 
class  RGB2XYZFunctor
 Convert linear (raw) RGB into standardized tri-stimulus XYZ. More...
 
class  RGBAccessor
 
class  RGBGradientMagnitudeFunctor
 
class  RGBPrime2LabFunctor
 Convert non-linear (gamma corrected) R'G'B' into perceptual uniform CIE L*a*b*. More...
 
class  RGBPrime2LuvFunctor
 Convert non-linear (gamma corrected) R'G'B' into perceptual uniform CIE L*u*v*. More...
 
class  RGBPrime2RGBFunctor
 Convert non-linear (gamma corrected) R'G'B' into non-linear (raw) RGB. More...
 
class  RGBPrime2XYZFunctor
 Convert non-linear (gamma corrected) R'G'B' into standardized tri-stimulus XYZ. More...
 
class  RGBPrime2YPrimeCbCrFunctor
 Convert non-linear (gamma corrected) R'G'B' into Y'CbCr color difference components. More...
 
class  RGBPrime2YPrimeIQFunctor
 Convert non-linear (gamma corrected) R'G'B' into Y'IQ components. More...
 
class  RGBPrime2YPrimePbPrFunctor
 Convert non-linear (gamma corrected) R'G'B' into Y'PbPr color difference components. More...
 
class  RGBPrime2YPrimeUVFunctor
 Convert non-linear (gamma corrected) R'G'B' into Y'UV components. More...
 
class  RGBToGrayAccessor
 
class  RGBValue
 Class for a single RGB value. More...
 
class  RowIterator
 Iterator adapter to linearly access row. More...
 
class  Sampler
 Create random samples from a sequence of indices. More...
 
class  SamplerOptions
 Options object for the Sampler class. More...
 
class  SeedOptions
 Options object for generateWatershedSeeds(). More...
 
class  SeedRgDirectValueFunctor
 Statistics functor to be used for seeded region growing. More...
 
class  SequenceAccessor
 Accessor for items that are STL compatible sequences. More...
 
class  ShortestPathDijkstra
 shortest path computer More...
 
class  SIFImportInfo
 Extracts image properties from an Andor SIF file header. More...
 
class  Size2D
 Two dimensional size object. More...
 
struct  SkeletonOptions
 Option object for skeletonizeImage() More...
 
class  SlantedEdgeMTFOptions
 Pass options to one of the slantedEdgeMTF() functions. More...
 
struct  SlicOptions
 Options object for slicSuperpixels(). More...
 
class  SortSamplesByDimensions
 
class  Splice
 
class  SplineImageView
 Create a continuous view onto a discrete image using splines. More...
 
class  SplineImageView0
 Create an image view for nearest-neighbor interpolation. More...
 
class  SplineImageView1
 Create an image view for bi-linear interpolation. More...
 
class  SplitBase
 
class  sRGB2RGBFunctor
 Convert standardized sRGB into non-linear (raw) RGB. More...
 
class  StandardAccessor
 Encapsulate access to the values an iterator points to. More...
 
class  StandardConstAccessor
 Encapsulate read access to the values an iterator points to. More...
 
class  StandardConstValueAccessor
 Encapsulate access to the values an iterator points to. More...
 
class  StandardValueAccessor
 Encapsulate access to the values an iterator points to. More...
 
class  StarPolytope
 Specialization of the polytope to polytopes which forms a star domain. More...
 
class  StaticPolynomial
 
class  StopAfterTree
 
class  StopAfterVoteCount
 
class  StopBase
 
class  StopIfBinTest
 
class  StopIfConverging
 
class  StopIfMargin
 
class  StopIfProb
 
struct  StridedArrayTag
 
class  StridedImageIterator
 Iterator to be used when pixels are to be skipped. More...
 
class  StridedMultiIterator
 A multi-dimensional hierarchical iterator to be used with vigra::MultiArrayView if it is not strided. More...
 
class  StridedScanOrderIterator
 Sequential iterator for MultiArrayView. More...
 
struct  ThinPlateSplineFunctor
 
class  ThreadPool
 Thread pool class to manage a set of parallel workers. More...
 
class  Threshold
 Threshold an image. More...
 
class  ThresholdSplit
 
class  TinyVector
 Class for fixed size vectors.This class contains an array of size SIZE of the specified VALUETYPE. The interface conforms to STL vector, except that there are no functions that change the size of a TinyVector. More...
 
class  TinyVectorBase
 Base class for fixed size vectors. More...
 
class  TinyVectorView
 Wrapper for fixed size vectors. More...
 
class  UniformIntRandomFunctor
 
class  UniformRandomFunctor
 
struct  UnstridedArrayTag
 
class  VectorAccessor
 Accessor for items that are STL compatible vectors. More...
 
class  VectorComponentAccessor
 Accessor for one component of a vector. More...
 
class  VectorComponentValueAccessor
 Accessor for one component of a vector. More...
 
class  VectorElementAccessor
 Accessor for one component of a vector. More...
 
class  VectorNormFunctor
 A functor for computing the vector norm. More...
 
class  VectorNormSqFunctor
 A functor for computing the squared vector norm. More...
 
class  VolumeExportInfo
 Argument object for the function exportVolume(). More...
 
class  VolumeImportInfo
 Argument object for the function importVolume(). More...
 
class  WatershedOptions
 Options object for watershed algorithms. More...
 
class  WignerMatrix
 computation of Wigner D matrix + rotation functions in SH,VH and R³ More...
 
class  XYZ2LabFunctor
 Convert standardized tri-stimulus XYZ into perceptual uniform CIE L*a*b*. More...
 
class  XYZ2LuvFunctor
 Convert standardized tri-stimulus XYZ into perceptual uniform CIE L*u*v*. More...
 
class  XYZ2RGBFunctor
 Convert standardized tri-stimulus XYZ into linear (raw) RGB. More...
 
class  XYZ2RGBPrimeFunctor
 Convert standardized tri-stimulus XYZ into non-linear (gamma corrected) R'G'B'. More...
 
class  YPrimeCbCr2RGBPrimeFunctor
 Convert Y'CbCr color difference components into non-linear (gamma corrected) R'G'B'. More...
 
class  YPrimeIQ2RGBPrimeFunctor
 Convert Y'IQ color components into non-linear (gamma corrected) R'G'B'. More...
 
class  YPrimePbPr2RGBPrimeFunctor
 Convert Y'PbPr color difference components into non-linear (gamma corrected) R'G'B'. More...
 
class  YPrimeUV2RGBPrimeFunctor
 Convert Y'UV color components into non-linear (gamma corrected) R'G'B'. More...
 

Typedefs

typedef AtImageBorder AtVolumeBorder
 Encode whether a voxel is near the volume border. More...
 
typedef BasicImage< UInt8BImage
 
typedef BasicImage< RGBValue
< UInt8 > > 
BRGBImage
 
typedef Diff2D CoordinateIterator
 Simulate an image where each pixel contains its coordinate. More...
 
typedef vigra::TinyVector< int, 3 > Diff3D
 3-dimensional difference vector
 
typedef BasicImage< double > DImage
 
typedef BasicImage< RGBValue
< double > > 
DRGBImage
 
typedef BasicImage< TinyVector
< double, 2 > > 
DVector2Image
 
typedef BasicImage< TinyVector
< double, 3 > > 
DVector3Image
 
typedef BasicImage< TinyVector
< double, 4 > > 
DVector4Image
 
typedef
EightNeighborhood::NeighborCode 
EightNeighborCode
 
typedef
NeighborOffsetCirculator
< EightNeighborCode
EightNeighborOffsetCirculator
 
typedef BasicImage< FFTWComplex<> > FFTWComplexImage
 
typedef BasicImage< fftw_real > FFTWRealImage
 
typedef BasicImage< float > FImage
 
typedef
FourNeighborhood::NeighborCode 
FourNeighborCode
 
typedef
NeighborOffsetCirculator
< FourNeighborCode
FourNeighborOffsetCirculator
 
typedef BasicImage< RGBValue
< float > > 
FRGBImage
 
typedef BasicImage< TinyVector
< float, 2 > > 
FVector2Image
 
typedef BasicImage< TinyVector
< float, 3 > > 
FVector3Image
 
typedef BasicImage< TinyVector
< float, 4 > > 
FVector4Image
 
typedef RidgeSplit
< BestGiniOfColumn
< GiniCriterion > > 
GiniRidgeSplit
 
typedef BasicImage< Int32IImage
 
typedef
detail::SelectIntegerType
< 16, detail::SignedIntTypes >
::type 
Int16
 16-bit signed int
 
typedef BasicImage< Int16Int16Image
 
typedef BasicImage< RGBValue
< Int16 > > 
Int16RGBImage
 
typedef
detail::SelectIntegerType
< 32, detail::SignedIntTypes >
::type 
Int32
 32-bit signed int
 
typedef BasicImage< Int32Int32Image
 
typedef BasicImage< RGBValue
< Int32 > > 
Int32RGBImage
 
typedef
detail::SelectIntegerType
< 64, detail::SignedIntTypes >
::type 
Int64
 64-bit signed int
 
typedef
detail::SelectIntegerType
< 8, detail::SignedIntTypes >
::type 
Int8
 8-bit signed int
 
typedef BasicImage< Int8Int8Image
 
typedef BasicImage< RGBValue
< Int8 > > 
Int8RGBImage
 
typedef
detail::SelectBiggestIntegerType
< detail::SignedIntTypes >
::type 
IntBiggest
 the biggest signed integer type of the system
 
typedef BasicImage< RGBValue
< Int32 > > 
IRGBImage
 
typedef RandomNumberGenerator
< detail::RandomState
< detail::MT19937 > > 
MersenneTwister
 
typedef std::ptrdiff_t MultiArrayIndex
 
typedef
Neighborhood3DSix::NeighborCode3D 
NeighborCode3DSix
 
typedef
Neighborhood3DTwentySix::NeighborCode3D 
NeighborCode3DTwentySix
 
typedef RandomNumberGenerator
< detail::RandomState
< detail::MT19937 > > 
RandomMT19937
 
typedef RandomNumberGenerator
< detail::RandomState
< detail::TT800 > > 
RandomTT800
 
typedef BasicImage< Int16SImage
 
typedef BasicImage< RGBValue
< Int16 > > 
SRGBImage
 
typedef RandomNumberGenerator
< detail::RandomState
< detail::TT800 > > 
TemperedTwister
 
typedef
detail::SelectIntegerType
< 16, detail::UnsignedIntTypes >
::type 
UInt16
 16-bit unsigned int
 
typedef BasicImage< UInt16UInt16Image
 
typedef BasicImage< RGBValue
< UInt16 > > 
UInt16RGBImage
 
typedef
detail::SelectIntegerType
< 32, detail::UnsignedIntTypes >
::type 
UInt32
 32-bit unsigned int
 
typedef BasicImage< UInt32UInt32Image
 
typedef BasicImage< RGBValue
< UInt32 > > 
UInt32RGBImage
 
typedef
detail::SelectIntegerType
< 64, detail::UnsignedIntTypes >
::type 
UInt64
 64-bit unsigned int
 
typedef
detail::SelectIntegerType
< 8, detail::UnsignedIntTypes >
::type 
UInt8
 8-bit unsigned int
 
typedef BasicImage< UInt8UInt8Image
 
typedef BasicImage< RGBValue
< UInt8 > > 
UInt8RGBImage
 
typedef
detail::SelectBiggestIntegerType
< detail::UnsignedIntTypes >
::type 
UIntBiggest
 the biggest unsigned integer type of the system
 

Enumerations

enum  AtImageBorder {
  NotAtBorder = 0, RightBorder = 1, LeftBorder = 2, TopBorder = 4,
  BottomBorder = 8, FrontBorder = 16
}
 Encode whether a point is near the image border. More...
 
enum  BoundaryDistanceTag { OuterBoundary, InterpixelBoundary, InnerBoundary }
 Specify which boundary is used for boundaryMultiDistance(). More...
 
enum  MultiArrayInitializationTag { LinearSequence }
 Initialize a MultiArray in a standard way. More...
 
enum  NeighborhoodType { DirectNeighborhood =0, IndirectNeighborhood =1 }
 Choose the neighborhood system in a dimension-independent way. More...
 
enum  RF_OptionTag
 
enum  SRGType
 

Functions

template<typename IntType >
Rational< IntType > abs (const Rational< IntType > &r)
 absolute value
 
template<typename Type >
Quaternion< Type >::NormType abs (Quaternion< Type > const &q)
 norm
 
template<unsigned IntBits, unsigned FracBits>
FixedPoint< IntBits, FracBits > abs (FixedPoint< IntBits, FracBits > v)
 absolute value.
 
template<class T , unsigned int RIDX, unsigned int GIDX, unsigned int BIDX>
RGBValue< T, RIDX, GIDX, BIDX > abs (RGBValue< T, RIDX, GIDX, BIDX > const &v)
 component-wise absolute value
 
template<class R >
FFTWComplex< R >::NormType abs (const FFTWComplex< R > &a)
 absolute value (= magnitude)
 
template<int IntBits, FPOverflowHandling OverflowHandling>
FixedPoint16< IntBits,
OverflowHandling > 
abs (FixedPoint16< IntBits, OverflowHandling > v)
 absolute value.
 
template<class V , int SIZE, class D1 , class D2 >
TinyVector< V, SIZE > abs (TinyVectorBase< V, SIZE, D1, D2 > const &v)
 component-wise absolute value
 
template<unsigned IntBits1, unsigned FracBits1, unsigned IntBits2, unsigned FracBits2, unsigned IntBits3, unsigned FracBits3>
void add (FixedPoint< IntBits1, FracBits1 > l, FixedPoint< IntBits2, FracBits2 > r, FixedPoint< IntBits3, FracBits3 > &result)
 addition with enforced result type.
 
template<int IntBits1, FPOverflowHandling OverflowHandling, int IntBits2, int IntBits3>
FixedPoint16< IntBits3,
OverflowHandling > & 
add (FixedPoint16< IntBits1, OverflowHandling > l, FixedPoint16< IntBits2, OverflowHandling > r, FixedPoint16< IntBits3, OverflowHandling > &result)
 addition with enforced result type.
 
template<class SrcIterator , class DestIterator >
linalg::TemporaryMatrix< double > affineMatrix2DFromCorrespondingPoints (SrcIterator s, SrcIterator send, DestIterator d)
 Create homogeneous matrix that maps corresponding points onto each other. More...
 
template<... >
void affineWarpImage (...)
 Warp an image according to an affine transformation. More...
 
template<class V1 , int SIZE, class D1 , class D2 , class V2 , class D3 , class D4 >
bool allGreater (TinyVectorBase< V1, SIZE, D1, D2 > const &l, TinyVectorBase< V2, SIZE, D3, D4 > const &r)
 pointwise greater-than
 
template<class V1 , int SIZE, class D1 , class D2 , class V2 , class D3 , class D4 >
bool allGreaterEqual (TinyVectorBase< V1, SIZE, D1, D2 > const &l, TinyVectorBase< V2, SIZE, D3, D4 > const &r)
 pointwise greater-equal
 
template<class V1 , int SIZE, class D1 , class D2 , class V2 , class D3 , class D4 >
bool allLess (TinyVectorBase< V1, SIZE, D1, D2 > const &l, TinyVectorBase< V2, SIZE, D3, D4 > const &r)
 pointwise less-than
 
template<class V1 , int SIZE, class D1 , class D2 , class V2 , class D3 , class D4 >
bool allLessEqual (TinyVectorBase< V1, SIZE, D1, D2 > const &l, TinyVectorBase< V2, SIZE, D3, D4 > const &r)
 pointwise less-equal
 
double angularGaborSigma (int directionCount, double centerFrequency)
 Calculate sensible angular sigma for given parameters. More...
 
template<... >
void anisotropicTotalVariationFilter (...)
 Performs Anisotropic Total Variation Regularization. More...
 
template<... >
void applyFourierFilter (...)
 Apply a filter (defined in the frequency domain) to an image. More...
 
template<... >
void applyFourierFilterFamily (...)
 Apply an array of filters (defined in the frequency domain) to an image. More...
 
template<class IndexIterator , class InIterator , class OutIterator >
void applyPermutation (IndexIterator index_first, IndexIterator index_last, InIterator in, OutIterator out)
 Sort an array according to the given index permutation. More...
 
template<class R >
arg (const FFTWComplex< R > &a)
 phase
 
template<class Iterator >
Iterator argMax (Iterator first, Iterator last)
 Find the maximum element in a sequence. More...
 
template<class Iterator , class UnaryFunctor >
Iterator argMaxIf (Iterator first, Iterator last, UnaryFunctor condition)
 Find the maximum element in a sequence conforming to a condition. More...
 
template<class Iterator >
Iterator argMin (Iterator first, Iterator last)
 Find the minimum element in a sequence. More...
 
template<class Iterator , class UnaryFunctor >
Iterator argMinIf (Iterator first, Iterator last, UnaryFunctor condition)
 Find the minimum element in a sequence conforming to a condition. More...
 
template<class T >
std::string asString (T t)(...)
 
template<int IntBits, FPOverflowHandling OverflowHandling>
FixedPoint16< 2, OverflowHandling > atan2 (FixedPoint16< IntBits, OverflowHandling > y, FixedPoint16< IntBits, OverflowHandling > x)
 Arctangent. Accuracy better than 1/3 degree (9 significant bits).
 
template<... >
void beaudetCornerDetector (...)
 Find corners in an image (4). More...
 
template<... >
void beautifyCrackEdgeImage (...)
 Beautify crack edge image for visualization. More...
 
double besselJ (int n, double x)
 Bessel function of the first kind. More...
 
double besselY (int n, double x)
 Bessel function of the second kind. More...
 
template<... >
void boundaryMultiDistance (...)
 Euclidean distance to the implicit boundaries of a multi-dimensional label array. More...
 
template<... >
void boundaryTensor (...)
 Calculate the boundary tensor for a scalar valued image. More...
 
template<... >
void boundaryTensor1 (...)
 Boundary tensor variant. More...
 
template<... >
void boundaryVectorDistance (...)
 Compute the vector distance transform to the implicit boundaries of a multi-dimensional label array. More...
 
template<... >
void cannyEdgeImage (...)
 Detect and mark edges in an edge image using Canny's algorithm. More...
 
template<... >
void cannyEdgeImageFromGradWithThinning (...)
 Detect and mark edges in an edge image using Canny's algorithm. More...
 
template<... >
void cannyEdgeImageWithThinning (...)
 Detect and mark edges in an edge image using Canny's algorithm. More...
 
template<... >
void cannyEdgelList (...)
 Simple implementation of Canny's edge detector. More...
 
template<... >
void cannyEdgelList3x3 (...)
 Improved implementation of Canny's edge detector. More...
 
template<... >
void cannyEdgelList3x3Threshold (...)
 Improved implementation of Canny's edge detector with thresholding. More...
 
template<... >
void cannyEdgelListThreshold (...)
 Canny's edge detector with thresholding. More...
 
template<class GRAPH , class EDGE_WEIGHTS , class SEEDS , class LABELS >
void carvingSegmentation (const GRAPH &g, const EDGE_WEIGHTS &edgeWeights, const SEEDS &seeds, const typename LABELS::Value backgroundLabel, const typename EDGE_WEIGHTS::Value backgroundBias, const typename EDGE_WEIGHTS::Value noPriorBelow, LABELS &labels)
 edge weighted watersheds Segmentataion More...
 
template<typename IntType >
Rational< IntType > ceil (const Rational< IntType > &r)
 smallest integer not smaller than r
 
template<unsigned IntBits, unsigned FracBits>
int ceil (FixedPoint< IntBits, FracBits > v)
 rounding up.
 
template<class V , unsigned int RIDX, unsigned int GIDX, unsigned int BIDX>
RGBValue< V, RIDX, GIDX, BIDX > ceil (RGBValue< V, RIDX, GIDX, BIDX > const &r)
 
template<int IntBits, FPOverflowHandling OverflowHandling>
Int32 ceil (FixedPoint16< IntBits, OverflowHandling > v)
 rounding up.
 
template<class V , int SIZE, class D1 , class D2 >
TinyVector< V, SIZE > ceil (TinyVectorBase< V, SIZE, D1, D2 > const &v)
 
UInt32 ceilPower2 (UInt32 x)
 Round up to the nearest power of 2. More...
 
template<int SIZE>
TinyVector< UInt32, SIZE > ceilPower2 (vigra::TinyVector< UInt32, SIZE > const &t)
 Round up values to the nearest power of 2. Implemented only for UInt32.
 
UInt32 checksum (const char *data, unsigned int size)
 Compute the CRC-32 checksum of a byte array. More...
 
double chi2 (unsigned int degreesOfFreedom, double arg, double accuracy=1e-7)
 Chi square distribution. More...
 
double chi2CDF (unsigned int degreesOfFreedom, double arg, double accuracy=1e-7)
 Cumulative chi square distribution. More...
 
template<class V , int SIZE>
TinyVector< V, SIZE > clip (TinyVector< V, SIZE > const &t, const V valLower, const V valUpper)
 Clip values to an interval. More...
 
template<class V , int SIZE>
TinyVector< V, SIZE > clip (TinyVector< V, SIZE > const &t, TinyVector< V, SIZE > const &valLower, TinyVector< V, SIZE > const &valUpper)
 Clip values to a vector of intervals. More...
 
template<class V , int SIZE>
TinyVector< V, SIZE > clipLower (TinyVector< V, SIZE > const &t)
 Clip negative values. More...
 
template<class V , int SIZE>
TinyVector< V, SIZE > clipLower (TinyVector< V, SIZE > const &t, const V val)
 Clip values below a threshold. More...
 
template<class V , int SIZE>
TinyVector< V, SIZE > clipUpper (TinyVector< V, SIZE > const &t, const V val)
 Clip values above a threshold. More...
 
template<class T1 , class T2 >
bool closeAtTolerance (T1 l, T2 r, typename PromoteTraits< T1, T2 >::Promote epsilon)
 Tolerance based floating-point equality. More...
 
template<... >
void closeGapsInCrackEdgeImage (...)
 Close one-pixel wide gaps in a cell grid edge image. More...
 
template<... >
void combineThreeImages (...)
 Combine three source images into destination image. More...
 
template<... >
void combineThreeMultiArrays (...)
 Combine three multi-dimensional arrays into one using a ternary function or functor. More...
 
template<... >
void combineTwoImages (...)
 Combine two source images into destination image. More...
 
template<... >
void combineTwoImagesIf (...)
 Combine ROI of two source images into destination image. More...
 
template<... >
void combineTwoMultiArrays (...)
 Combine two multi-dimensional arrays into one using a binary function or functor. More...
 
void compress (char const *source, std::size_t size, ArrayVector< char > &dest, CompressionMethod method)
 
UInt32 concatenateChecksum (UInt32 checksum, const char *data, unsigned int size)
 
template<class ValueType >
Quaternion< ValueType > conj (Quaternion< ValueType > const &q)
 Create conjugate quaternion.
 
template<class R >
FFTWComplex< R > conj (const FFTWComplex< R > &a)
 complex conjugate
 
template<class PointArray1 , class PointArray2 >
void convexHull (const PointArray1 &points, PointArray2 &convex_hull)
 Compute convex hull of a 2D polygon. More...
 
template<... >
void convolveFFT (...)
 Convolve an array with a kernel by means of the Fourier transform. More...
 
template<... >
void convolveFFTComplex (...)
 Convolve a complex-valued array by means of the Fourier transform. More...
 
template<... >
void convolveFFTComplexMany (...)
 Convolve a complex-valued array with a sequence of kernels by means of the Fourier transform. More...
 
template<... >
void convolveFFTMany (...)
 Convolve a real-valued array with a sequence of kernels by means of the Fourier transform. More...
 
template<... >
void convolveImage (...)
 Convolve an image with the given kernel(s). More...
 
template<... >
void convolveImageWithMask (...)
 Deprecated name of 2-dimensional normalized convolution, i.e. convolution with a mask image. More...
 
template<... >
void convolveLine (...)
 Performs a 1-dimensional convolution of the source signal using the given kernel. More...
 
template<... >
void convolveMultiArrayOneDimension (...)
 Convolution along a single dimension of a multi-dimensional arrays. More...
 
template<class G , class A , class B >
void copyEdgeMap (const G &g, const A &a, B &b)
 copy a lemon edge map
 
template<... >
void copyImage (...)
 Copy source image into destination image. More...
 
template<... >
void copyImageIf (...)
 Copy source ROI into destination image. More...
 
template<... >
void copyMultiArray (...)
 Copy a multi-dimensional array. More...
 
template<class G , class A , class B >
void copyNodeMap (const G &g, const A &a, B &b)
 copy a lemon node map
 
template<... >
void cornerResponseFunction (...)
 Find corners in an image (1). More...
 
template<... >
void correlateFFT (...)
 Correlate an array with a kernel by means of the Fourier transform. More...
 
template<class REAL >
REAL cos_pi (REAL x)
 cos(pi*x). More...
 
template<int N>
CoupledIteratorType< N >::type createCoupledIterator (TinyVector< MultiArrayIndex, N > const &shape)
 
template<unsigned int N1, class T1 , class S1 >
CoupledIteratorType< N1, T1 >::type createCoupledIterator (MultiArrayView< N1, T1, S1 > const &m1)
 
template<unsigned int N1, class T1 , class S1 , unsigned int N2, class T2 , class S2 >
CoupledIteratorType< N1, T1,
T2 >::type 
createCoupledIterator (MultiArrayView< N1, T1, S1 > const &m1, MultiArrayView< N2, T2, S2 > const &m2)
 
template<unsigned int N1, class T1 , class S1 , unsigned int N2, class T2 , class S2 , unsigned int N3, class T3 , class S3 >
CoupledIteratorType< N1, T1,
T2, T3 >::type 
createCoupledIterator (MultiArrayView< N1, T1, S1 > const &m1, MultiArrayView< N2, T2, S2 > const &m2, MultiArrayView< N3, T3, S3 > const &m3)
 
template<unsigned int N1, class T1 , class S1 , unsigned int N2, class T2 , class S2 , unsigned int N3, class T3 , class S3 , unsigned int N4, class T4 , class S4 >
CoupledIteratorType< N1, T1,
T2, T3, T4 >::type 
createCoupledIterator (MultiArrayView< N1, T1, S1 > const &m1, MultiArrayView< N2, T2, S2 > const &m2, MultiArrayView< N3, T3, S3 > const &m3, MultiArrayView< N4, T4, S4 > const &m4)
 
template<unsigned int N1, class T1 , class S1 , unsigned int N2, class T2 , class S2 , unsigned int N3, class T3 , class S3 , unsigned int N4, class T4 , class S4 , unsigned int N5, class T5 , class S5 >
CoupledIteratorType< N1, T1,
T2, T3, T4, T5 >::type 
createCoupledIterator (MultiArrayView< N1, T1, S1 > const &m1, MultiArrayView< N2, T2, S2 > const &m2, MultiArrayView< N3, T3, S3 > const &m3, MultiArrayView< N4, T4, S4 > const &m4, MultiArrayView< N5, T5, S5 > const &m5)
 
template<unsigned int N, class T >
ChunkedArray< N, T >::iterator createCoupledIterator (ChunkedArray< N, T > &m)
 
template<... >
void createGaborFilter (...)
 Create a gabor filter in frequency space. More...
 
template<int N>
HierarchicalIteratorType< N >::type createHierarchicalIterator (TinyVector< MultiArrayIndex, N > const &shape)
 
template<unsigned int N1, class T1 , class S1 >
HierarchicalIteratorType< N1,
T1 >::type 
createHierarchicalIterator (MultiArrayView< N1, T1, S1 > const &m1)
 
template<unsigned int N1, class T1 , class S1 , unsigned int N2, class T2 , class S2 >
HierarchicalIteratorType< N1,
T1, T2 >::type 
createHierarchicalIterator (MultiArrayView< N1, T1, S1 > const &m1, MultiArrayView< N2, T2, S2 > const &m2)
 
template<unsigned int N1, class T1 , class S1 , unsigned int N2, class T2 , class S2 , unsigned int N3, class T3 , class S3 >
HierarchicalIteratorType< N1,
T1, T2, T3 >::type 
createHierarchicalIterator (MultiArrayView< N1, T1, S1 > const &m1, MultiArrayView< N2, T2, S2 > const &m2, MultiArrayView< N3, T3, S3 > const &m3)
 
template<unsigned int N1, class T1 , class S1 , unsigned int N2, class T2 , class S2 , unsigned int N3, class T3 , class S3 , unsigned int N4, class T4 , class S4 >
HierarchicalIteratorType< N1,
T1, T2, T3, T4 >::type 
createHierarchicalIterator (MultiArrayView< N1, T1, S1 > const &m1, MultiArrayView< N2, T2, S2 > const &m2, MultiArrayView< N3, T3, S3 > const &m3, MultiArrayView< N4, T4, S4 > const &m4)
 
template<unsigned int N1, class T1 , class S1 , unsigned int N2, class T2 , class S2 , unsigned int N3, class T3 , class S3 , unsigned int N4, class T4 , class S4 , unsigned int N5, class T5 , class S5 >
HierarchicalIteratorType< N1,
T1, T2, T3, T4, T5 >::type 
createHierarchicalIterator (MultiArrayView< N1, T1, S1 > const &m1, MultiArrayView< N2, T2, S2 > const &m2, MultiArrayView< N3, T3, S3 > const &m3, MultiArrayView< N4, T4, S4 > const &m4, MultiArrayView< N5, T5, S5 > const &m5)
 
template<... >
void createRGBTiffImage (...)
 Create a 3-band TiffImage from the given RGB image. More...
 
template<... >
void createScalarTiffImage (...)
 Create a single-band TiffImage from the given scalar image. More...
 
template<... >
void createTiffImage (...)
 Create a TiffImage from the given iterator range. More...
 
template<class V1 , unsigned int R, unsigned int G, unsigned int B, class V2 >
PromoteTraits< RGBValue< V1, R,
G, B >, RGBValue< V2, R, G, B >
>::Promote 
cross (RGBValue< V1, R, G, B > const &r1, RGBValue< V2, R, G, B > const &r2)
 cross product
 
template<class V1 , class D1 , class D2 , class V2 , class D3 , class D4 >
TinyVector< typename
PromoteTraits< V1, V2 >
::Promote, 3 > 
cross (TinyVectorBase< V1, 3, D1, D2 > const &r1, TinyVectorBase< V2, 3, D3, D4 > const &r2)
 cross product
 
template<... >
void crossCorrelation (...)
 This function performes a (slow) cross-correlation. More...
 
template<class V , int SIZE, class D1 , class D2 >
TinyVector< typename
NumericTraits< V >::Promote,
SIZE > 
cumprod (TinyVectorBase< V, SIZE, D1, D2 > const &l)
 cumulative product of the vector's elements
 
template<class V , int SIZE, class D1 , class D2 >
TinyVector< typename
NumericTraits< V >::Promote,
SIZE > 
cumsum (TinyVectorBase< V, SIZE, D1, D2 > const &l)
 cumulative sum of the vector's elements
 
template<... >
void differenceOfExponentialCrackEdgeImage (...)
 Detect and mark edges in a crack edge image using the Shen/Castan zero-crossing detector. More...
 
template<... >
void differenceOfExponentialEdgeImage (...)
 Detect and mark edges in an edge image using the Shen/Castan zero-crossing detector. More...
 
template<... >
void discDilation (...)
 Apply dilation (maximum) filter with disc of given radius to image. More...
 
template<... >
void discDilationWithMask (...)
 Apply dilation (maximum) filter with disc of given radius to image using a mask. More...
 
template<... >
void discErosion (...)
 Apply erosion (minimum) filter with disc of given radius to image. More...
 
template<... >
void discErosionWithMask (...)
 Apply erosion (minimum) filter with disc of given radius to image using a mask. More...
 
template<... >
void discMedian (...)
 Apply median filter with disc of given radius to image. More...
 
template<... >
void discMedianWithMask (...)
 Apply median filter with disc of given radius to image using a mask. More...
 
template<... >
void discRankOrderFilter (...)
 Apply rank order filter with disc structuring function to the image. More...
 
template<... >
void discRankOrderFilterWithMask (...)
 Apply rank order filter with disc structuring function to the image using a mask. More...
 
template<... >
void distanceTransform (...)
 
template<int IntBits1, FPOverflowHandling OverflowHandling, int IntBits2, int IntBits3>
FixedPoint16< IntBits3,
OverflowHandling > & 
div (FixedPoint16< IntBits1, OverflowHandling > l, FixedPoint16< IntBits2, OverflowHandling > r, FixedPoint16< IntBits3, OverflowHandling > &result)
 division with enforced result type.
 
template<class V , int SIZE, class D1 , class D2 >
TinyVector< V, SIZE > div (TinyVectorBase< V, SIZE, D1, D2 > const &l, V v)
 component-wise scalar division without type promotion
 
template<class V1 , unsigned int RIDX1, unsigned int GIDX1, unsigned int BIDX1, class V2 , unsigned int RIDX2, unsigned int GIDX2, unsigned int BIDX2>
PromoteTraits< V1, V2 >::Promote dot (RGBValue< V1, RIDX1, GIDX1, BIDX1 > const &r1, RGBValue< V2, RIDX2, GIDX2, BIDX2 > const &r2)
 dot product
 
template<class V1 , int SIZE, class D1 , class D2 , class V2 , class D3 , class D4 >
PromoteTraits< V1, V2 >::Promote dot (TinyVectorBase< V1, SIZE, D1, D2 > const &l, TinyVectorBase< V2, SIZE, D3, D4 > const &r)
 dot product
 
 doxygen_overloaded_function (template<...> void separableConvolveBlockwise) template< unsigned int N
 Separated convolution on ChunkedArrays. More...
 
template<unsigned IntBits, unsigned FracBits>
FixedPoint< 0, FracBits > dual_frac (FixedPoint< IntBits, FracBits > v)
 dual fractional part: 1 - frac(v).
 
template<int IntBits, FPOverflowHandling OverflowHandling>
FixedPoint16< IntBits,
OverflowHandling > 
dual_frac (FixedPoint16< IntBits, OverflowHandling > v)
 dual fractional part. (1 - frac(v))
 
template<unsigned int N, class T , class S , class Array >
void eccentricityCenters (const MultiArrayView< N, T, S > &src, Array &centers)
 Find the (approximate) eccentricity center in each region of a labeled image. More...
 
template<unsigned int N, class T , class S , class Array >
void eccentricityTransformOnLabels (MultiArrayView< N, T > const &src, MultiArrayView< N, S > dest, Array &centers)
 Computes the (approximate) eccentricity transform on each region of a labeled image. More...
 
template<class GRAPH , class WEIGHTS , class COMPERATOR >
void edgeSort (const GRAPH &g, const WEIGHTS &weights, const COMPERATOR &comperator, std::vector< typename GRAPH::Edge > &sortedEdges)
 get a vector of Edge descriptors More...
 
template<class GRAPH , class EDGE_WEIGHTS , class SEEDS , class LABELS >
void edgeWeightedWatershedsSegmentation (const GRAPH &g, const EDGE_WEIGHTS &edgeWeights, const SEEDS &seeds, LABELS &labels)
 edge weighted watersheds Segmentataion More...
 
template<unsigned int N, class DirectedTag , class T , class EDGEMAP >
void edgeWeightsFromInterpolatedImage (const GridGraph< N, DirectedTag > &g, const MultiArrayView< N, T > &interpolatedImage, EDGEMAP &edgeWeights, bool euclidean=false)
 create edge weights from an interpolated image More...
 
template<unsigned int N, class DirectedTag , class NODEMAP , class EDGEMAP , class FUNCTOR >
void edgeWeightsFromNodeWeights (const GridGraph< N, DirectedTag > &g, const NODEMAP &nodeWeights, EDGEMAP &edgeWeights, bool euclidean, FUNCTOR const &func)
 create edge weights from node weights More...
 
double ellipticIntegralE (double x, double k)
 The incomplete elliptic integral of the second kind. More...
 
double ellipticIntegralF (double x, double k)
 The incomplete elliptic integral of the first kind. More...
 
template<... >
void estimateAffineTransform (...)
 Estimate the optical flow between two images according to an affine transform model (e.g. translation, rotation, non-uniform scaling, and shearing). More...
 
template<... >
void estimateGlobalRotation (...)
 Estimate the rotation between two images by means of a normalized cross correlation matching of the FFT spectra. More...
 
template<... >
void estimateGlobalRotationTranslation (...)
 Estimate the (global) rotation and translation between two images by means a normalized cross correlation matching. More...
 
template<... >
void estimateGlobalTranslation (...)
 Estimate the translation between two images by means of a normalized cross correlation matching. More...
 
template<... >
void estimateSimilarityTransform (...)
 Estimate the optical flow between two images according to a similarity transform model (e.g. translation, rotation, and uniform scaling). More...
 
template<... >
void estimateTranslation (...)
 Estimate the optical flow between two images according to a translation model. More...
 
bool even (int t)
 Check if an integer is even. More...
 
template<... >
void exportImage (...)
 Write an image to a file. More...
 
template<... >
void exportImageAlpha (...)
 Write the image and its alpha channel to a file. More...
 
template<... >
void exportVolume (...)
 Function for exporting a 3D volume. More...
 
template<... >
void extendedLocalMaxima (...)
 Find local maximal regions in an array. More...
 
template<... >
void extendedLocalMaxima3D (...)
 Find local maximal regions in 3D multi array. More...
 
template<... >
void extendedLocalMinima (...)
 Find local minimal regions (plateaus) in an array. More...
 
template<... >
void extendedLocalMinima3D (...)
 Find local minimal regions in a volume. More...
 
template<class T , class S , class PointArray >
void extractContour (MultiArrayView< 2, T, S > const &label_image, Shape2 const &anchor_point, PointArray &contour_points)
 Create a polygon from the interpixel contour of a labeled region. More...
 
template<... >
void fastCrossCorrelation (...)
 This function performes a fast cross-correlation. More...
 
template<... >
void fastNormalizedCrossCorrelation (...)
 This function performes a fast normalized cross-correlation. More...
 
template<class GRAPH , class EDGE_WEIGHTS , class NODE_SIZE , class NODE_LABEL_MAP >
void felzenszwalbSegmentation (const GRAPH &graph, const EDGE_WEIGHTS &edgeWeights, const NODE_SIZE &nodeSizes, float k, NODE_LABEL_MAP &nodeLabeling, const int nodeNumStopCond=-1)
 edge weighted watersheds Segmentataion More...
 
template<class T , int N>
TinyVector< T, N > fftwCorrespondingShapeR2C (TinyVector< T, N > shape)
 Find frequency domain shape for a R2C Fourier transform. More...
 
template<class G , class A , class T >
void fillEdgeMap (const G &g, A &a, const T &value)
 fill a lemon edge map
 
template<class G , class A , class T >
void fillNodeMap (const G &g, A &a, const T &value)
 fill a lemon node map
 
template<class Point , class T , class S , class Value >
void fillPolygon (Polygon< Point > const &p, MultiArrayView< 2, T, S > &output_image, Value value)
 Render closed polygon p into the image output_image. More...
 
template<class TARGET , unsigned IntBits, unsigned FracBits>
TARGET fixed_point_cast (FixedPoint< IntBits, FracBits > v)
 
template<class TARGET , int IntBits, FPOverflowHandling OverflowHandling>
TARGET fixed_point_cast (FixedPoint16< IntBits, OverflowHandling > v)
 
template<typename IntType >
Rational< IntType > floor (const Rational< IntType > &r)
 largest integer not larger than r
 
template<unsigned IntBits, unsigned FracBits>
int floor (FixedPoint< IntBits, FracBits > v)
 rounding down.
 
template<class V , unsigned int RIDX, unsigned int GIDX, unsigned int BIDX>
RGBValue< V, RIDX, GIDX, BIDX > floor (RGBValue< V, RIDX, GIDX, BIDX > const &r)
 
template<int IntBits, FPOverflowHandling OverflowHandling>
Int32 floor (FixedPoint16< IntBits, OverflowHandling > v)
 rounding down.
 
template<class V , int SIZE, class D1 , class D2 >
TinyVector< V, SIZE > floor (TinyVectorBase< V, SIZE, D1, D2 > const &v)
 
UInt32 floorPower2 (UInt32 x)
 Round down to the nearest power of 2. More...
 
template<int SIZE>
TinyVector< UInt32, SIZE > floorPower2 (vigra::TinyVector< UInt32, SIZE > const &t)
 Round down values to the nearest power of 2. Implemented only for UInt32.
 
template<... >
void foerstnerCornerDetector (...)
 Find corners in an image (2). More...
 
template<typename TPL , typename FUNCTOR >
void for_each_in_tuple (TPL &&t, FUNCTOR &&f)
 
template<... >
void fourierTransform (...)
 Compute forward and inverse Fourier transforms. More...
 
template<... >
void fourierTransformInverse (...)
 Compute inverse Fourier transforms. More...
 
template<... >
void fourierTransformReal (...)
 Real Fourier transforms for even and odd boundary conditions (aka. cosine and sine transforms). More...
 
template<unsigned IntBits, unsigned FracBits>
FixedPoint< 0, FracBits > frac (FixedPoint< IntBits, FracBits > v)
 fractional part.
 
template<int IntBits, FPOverflowHandling OverflowHandling>
FixedPoint16< IntBits,
OverflowHandling > 
frac (FixedPoint16< IntBits, OverflowHandling > v)
 fractional part. (difference between v and its floor)
 
double gamma (double x)
 The gamma function. More...
 
template<... >
void gaussianDivergenceMultiArray (...)
 Calculate the divergence of a vector field using Gaussian derivative filters. More...
 
template<... >
void gaussianGradient (...)
 Calculate the gradient vector by means of a 1st derivatives of Gaussian filter. More...
 
template<... >
void gaussianGradientMagnitude (...)
 Calculate the gradient magnitude by means of a 1st derivatives of Gaussian filter. More...
 
template<... >
void gaussianGradientMultiArray (...)
 Calculate Gaussian gradient of a multi-dimensional arrays. More...
 
template<... >
void gaussianSharpening (...)
 Perform sharpening function with gaussian filter. More...
 
template<... >
void gaussianSmoothing (...)
 Perform isotropic Gaussian convolution. More...
 
template<... >
void gaussianSmoothMultiArray (...)
 Isotropic Gaussian smoothing of a multi-dimensional arrays. More...
 
template<typename IntType >
IntType gcd (IntType n, IntType m)
 
template<... >
unsigned int generateSlicSeeds (...)
 Generate seeds for SLIC superpixel computation in arbitrary dimensions. More...
 
template<... >
unsigned int generateWatershedSeeds (...)
 Generate seeds for watershed computation and seeded region growing. More...
 
template<unsigned int TARGET_INDEX, class Handle >
CoupledHandleCast
< TARGET_INDEX, Handle >
::reference 
get (Handle &handle)
 
template<unsigned int TARGET_INDEX, class Handle >
CoupledHandleCast
< TARGET_INDEX, Handle >
::const_reference 
get (Handle const &handle)
 
template<... >
void getAnisotropy (...)
 Sets up directional data for anisotropic regularization. More...
 
template<... >
void gradientBasedTransform (...)
 Calculate a function of the image gradient. More...
 
template<... >
void gradientEnergyTensor (...)
 Calculate the gradient energy tensor for a scalar valued image. More...
 
template<class GRAPH , class NODE_FEATURES_IN , class EDGE_INDICATOR , class NODE_FEATURES_OUT >
void graphSmoothing (const GRAPH &g, const NODE_FEATURES_IN &nodeFeaturesIn, const EDGE_INDICATOR &edgeIndicator, const float lambda, const float edgeThreshold, const float scale, NODE_FEATURES_OUT &nodeFeaturesOut)
 smooth node features of a graph More...
 
template<class T1 , class T2 >
bool greaterEqualAtTolerance (T1 l, T2 r, typename PromoteTraits< T1, T2 >::Promote epsilon)
 Tolerance based floating-point greater-or-equal. More...
 
template<... >
void hessianMatrixOfGaussian (...)
 Filter image with the 2nd derivatives of the Gaussian at the given scale to get the Hessian matrix. More...
 
template<... >
void hessianOfGaussianMultiArray (...)
 Calculate Hessian matrix of a N-dimensional arrays using Gaussian derivative filters. More...
 
template<... >
void hierarchicalClustering (...)
 Reduce the number of nodes in a graph by iteratively contracting the cheapest edge. More...
 
template<... >
void hourGlassFilter (...)
 Anisotropic tensor smoothing with the hourglass filter. More...
 
template<int IntBits, FPOverflowHandling OverflowHandling>
FixedPoint16< IntBits,
OverflowHandling > 
hypot (FixedPoint16< IntBits, OverflowHandling > v1, FixedPoint16< IntBits, OverflowHandling > v2)
 Length of hypotenuse.
 
template<class R >
imag (const FFTWComplex< R > &a)
 imaginary part
 
std::string impexListExtensions ()
 List the file extension VIGRA understands. More...
 
std::string impexListFormats ()
 List the image formats VIGRA can read and write. More...
 
template<... >
void importImage (...)
 Read an image from a file. More...
 
template<... >
void importImageAlpha (...)
 Read the image specified by the given vigra::ImageImportInfo object including its alpha channel. More...
 
template<... >
void importTiffImage (...)
 Read a given TIFF image. More...
 
template<... >
void importVolume (...)
 Function for importing a 3D volume. More...
 
template<class Iterator , class IndexIterator , class Compare >
void indexSort (Iterator first, Iterator last, IndexIterator index_first, Compare c)
 Return the index permutation that would sort the input array. More...
 
template<... >
void initImage (...)
 Write a value to every pixel in an image or rectangular ROI. More...
 
template<... >
void initImageBorder (...)
 Write value to the specified border pixels in the image. More...
 
template<... >
void initImageIf (...)
 Write value to pixel in the image if mask is true. More...
 
template<... >
void initImageWithFunctor (...)
 Write the result of a functor call to every pixel in an image or rectangular ROI. More...
 
template<... >
void initMultiArray (...)
 Write a value to every element in a multi-dimensional array. More...
 
template<... >
void initMultiArrayBorder (...)
 Write values to the specified border values in the array. More...
 
template<... >
void inspectImage (...)
 Apply read-only functor to every pixel in the image. More...
 
template<... >
void inspectImageIf (...)
 Apply read-only functor to every pixel in the ROI. More...
 
template<... >
void inspectMultiArray (...)
 Call an analyzing functor at every element of a multi-dimensional array. More...
 
template<... >
void inspectSequence (...)
 Call an analyzing functor at every element of a sequence. More...
 
template<... >
void inspectTwoImages (...)
 Apply read-only functor to every pixel of both images. More...
 
template<... >
void inspectTwoImagesIf (...)
 Apply read-only functor to those pixels of both images where the mask image is non-zero. More...
 
template<... >
void inspectTwoMultiArrays (...)
 Call an analyzing functor at all corresponding elements of two multi-dimensional arrays. More...
 
template<class InIterator , class OutIterator >
void inversePermutation (InIterator first, InIterator last, OutIterator out)
 Compute the inverse of a given permutation. More...
 
AtImageBorder isAtImageBorder (int x, int y, int width, int height)
 Find out whether a point is at the image border. More...
 
AtVolumeBorder isAtVolumeBorder (int x, int y, int z, int width, int height, int depth)
 Find out whether a voxel is at the volume border. More...
 
AtVolumeBorder isAtVolumeBorderAntiCausal (int x, int y, int z, int width, int height, int depth)
 Find out whether a voxel is at a scan-order relevant volume border. This function checks if x == 0 or y == 0 or z == 0 and returns the appropriate value of vigra::AtVolumeBorder, or zero when the voxel is not at te volume border. The behavior of the function is undefined if (x,y,z) is not inside the volume. More...
 
AtVolumeBorder isAtVolumeBorderCausal (int x, int y, int z, int width, int height, int)
 Find out whether a voxel is at a scan-order relevant volume border. This function checks if x == 0 or y == 0 or z == 0 and returns the appropriate value of vigra::AtVolumeBorder, or zero when the voxel is not at te volume border. The behavior of the function is undefined if (x,y,z) is not inside the volume.
 
bool isHDF5 (char const *filename)
 Check if given filename refers to a HDF5 file.
 
bool isImage (char const *filename)
 Test whether a file is an image format known to VIGRA. More...
 
bool isPower2 (UInt32 x)
 Determine whether x is a power of 2 Bit twiddle from https://graphics.stanford.edu/~seander/bithacks.html#DetermineIfPowerOf2.
 
template<class V >
TinyVector< float, 3 > lab2Polar (V const &lab)
 Create polar representation form L*a*b*. More...
 
template<... >
unsigned int labelImage (...)
 Find the connected components of a segmented image. More...
 
template<... >
unsigned int labelImageWithBackground (...)
 Find the connected components of a segmented image, excluding the background from labeling. More...
 
template<... >
unsigned int labelMultiArray (...)
 Find the connected components of a MultiArray with arbitrary many dimensions. More...
 
template<... >
unsigned int labelMultiArrayBlockwise (...)
 Connected components labeling for MultiArrays and ChunkedArrays. More...
 
template<... >
unsigned int labelMultiArrayWithBackground (...)
 Find the connected components of a MultiArray with arbitrary many dimensions, excluding the background from labeling. More...
 
template<... >
unsigned int labelVolume (...)
 Find the connected components of a segmented volume. More...
 
template<class SrcIterator , class SrcAccessor , class SrcShape , class DestIterator , class DestAccessor >
unsigned int labelVolumeSix (triple< SrcIterator, SrcShape, SrcAccessor > src, pair< DestIterator, DestAccessor > dest)
 Find the connected components of a segmented volume using the 6-neighborhood. More...
 
template<... >
unsigned int labelVolumeWithBackground (...)
 Find the connected components of a segmented volume, excluding the background from labeling. More...
 
template<... >
void laplacianOfGaussian (...)
 Filter image with the Laplacian of Gaussian operator at the given scale. More...
 
template<... >
void laplacianOfGaussianMultiArray (...)
 Calculate Laplacian of a N-dimensional arrays using Gaussian derivative filters. More...
 
template<typename IntType >
IntType lcm (IntType n, IntType m)
 
template<class REAL >
REAL legendre (unsigned int l, int m, REAL x)
 Associated Legendre polynomial. More...
 
template<class REAL >
REAL legendre (unsigned int l, REAL x)
 Legendre polynomial. More...
 
template<class T1 , class T2 >
bool lessEqualAtTolerance (T1 l, T2 r, typename PromoteTraits< T1, T2 >::Promote epsilon)
 Tolerance based floating-point less-or-equal. More...
 
template<class Multiplier , class DestValueType >
LinearIntensityTransform
< DestValueType, Multiplier > 
linearIntensityTransform (Multiplier scale, DestValueType offset)
 Apply a linear transform to the source pixel values. More...
 
template<... >
bool linearNoiseNormalization (...)
 Noise normalization by means of an estimated or given linear noise model. More...
 
template<class SrcValueType , class DestValueType >
LinearIntensityTransform
< DestValueType, typename
NumericTraits< DestValueType >
::RealPromote > 
linearRangeMapping (SrcValueType src_min, SrcValueType src_max, DestValueType dest_min, DestValueType dest_max)
 Map a source intensity range linearly to a destination range. More...
 
template<class Iterator , class Value >
void linearSequence (Iterator first, Iterator last, Value start, Value step)
 Fill an array with a sequence of numbers. More...
 
template<... >
void localMaxima (...)
 Find local maxima in an image or multi-dimensional array. More...
 
template<... >
void localMaxima3D (...)
 Find local maxima in a 3D multi array. More...
 
template<... >
void localMinima (...)
 Find local minima in an image or multi-dimensional array. More...
 
template<... >
void localMinima3D (...)
 Find local minima in a 3D multi array. More...
 
Int32 log2i (UInt32 x)
 Compute the base-2 logarithm of an integer. More...
 
double loggamma (double x)
 The natural logarithm of the gamma function. More...
 
template<class V >
TinyVector< float, 3 > luv2Polar (V const &luv)
 Create polar representation form L*u*v*. More...
 
template<class T , class Stride >
BasicImageView< T > makeBasicImageView (MultiArrayView< 2, T, Stride > const &array)
 
template<class T >
BasicImageView< T > makeBasicImageView (MultiArray< 3, T > const &array)
 
template<class ArrayLike , class Compare >
detail::IndexCompare
< ArrayLike, Compare > 
makeIndexComparator (ArrayLike a, Compare c)
 Create a compare functor for indirect sort. More...
 
template<class GRAPH_IN , class GRAPH_IN_NODE_LABEL_MAP >
void makeRegionAdjacencyGraph (GRAPH_IN graphIn, GRAPH_IN_NODE_LABEL_MAP labels, AdjacencyListGraph &rag, typename AdjacencyListGraph::template EdgeMap< std::vector< typename GRAPH_IN::Edge > > &affiliatedEdges, const Int64 ignoreLabel=-1)
 make a region adjacency graph from a graph and labels w.r.t. that graph More...
 
template<class T , class Stride >
BasicImageView< RGBValue< T > > makeRGBImageView (MultiArrayView< 3, T, Stride > const &array)
 
template<class V1 , int SIZE, class D1 , class D2 , class V2 , class D3 , class D4 >
TinyVector< typename
PromoteTraits< V1, V2 >
::Promote, SIZE > 
max (TinyVectorBase< V1, SIZE, D1, D2 > const &l, TinyVectorBase< V2, SIZE, D3, D4 > const &r)
 element-wise maximum
 
template<class V , int SIZE, class D1 , class D2 >
V const & max (TinyVectorBase< V, SIZE, D1, D2 > const &l)
 maximum element
 
triple< Diff2D, Diff2D,
MeshGridAccessor
meshGrid (Diff2D upperLeft, Diff2D lowerRight)
 
template<class V1 , int SIZE, class D1 , class D2 , class V2 , class D3 , class D4 >
TinyVector< typename
PromoteTraits< V1, V2 >
::Promote, SIZE > 
min (TinyVectorBase< V1, SIZE, D1, D2 > const &l, TinyVectorBase< V2, SIZE, D3, D4 > const &r)
 element-wise minimum
 
template<class V , int SIZE, class D1 , class D2 >
V const & min (TinyVectorBase< V, SIZE, D1, D2 > const &l)
 minimum element
 
template<... >
void moveDCToCenter (...)
 Rearrange the quadrants of a Fourier image so that the origin is in the image center. More...
 
template<... >
void moveDCToUpperLeft (...)
 Rearrange the quadrants of a Fourier image so that the origin is in the image's upper left. More...
 
template<class Vector >
double mtfFitGaussian (Vector const &mtf)
 Fit a Gaussian function to a given MTF. More...
 
template<unsigned IntBits1, unsigned FracBits1, unsigned IntBits2, unsigned FracBits2, unsigned IntBits3, unsigned FracBits3>
void mul (FixedPoint< IntBits1, FracBits1 > l, FixedPoint< IntBits2, FracBits2 > r, FixedPoint< IntBits3, FracBits3 > &result)
 multiplication with enforced result type.
 
template<int IntBits1, FPOverflowHandling OverflowHandling, int IntBits2, int IntBits3>
FixedPoint16< IntBits3,
OverflowHandling > & 
mul (FixedPoint16< IntBits1, OverflowHandling > l, FixedPoint16< IntBits2, OverflowHandling > r, FixedPoint16< IntBits3, OverflowHandling > &result)
 multiplication with enforced result type.
 
template<... >
void multiBinaryDilation (...)
 Binary dilation on multi-dimensional arrays. More...
 
template<... >
void multiBinaryErosion (...)
 Binary erosion on multi-dimensional arrays. More...
 
template<... >
void multiGrayscaleDilation (...)
 Parabolic grayscale dilation on multi-dimensional arrays. More...
 
template<... >
void multiGrayscaleErosion (...)
 Parabolic grayscale erosion on multi-dimensional arrays. More...
 
template<... >
void noiseVarianceClustering (...)
 Determine the noise variance as a function of the image intensity and cluster the results. More...
 
template<... >
void noiseVarianceEstimation (...)
 Determine the noise variance as a function of the image intensity. More...
 
double noncentralChi2 (unsigned int degreesOfFreedom, double noncentrality, double arg, double accuracy=1e-7)
 Non-central chi square distribution. More...
 
double noncentralChi2CDF (unsigned int degreesOfFreedom, double noncentrality, double arg, double accuracy=1e-7)
 Cumulative non-central chi square distribution. More...
 
double noncentralChi2CDFApprox (unsigned int degreesOfFreedom, double noncentrality, double arg)
 Cumulative non-central chi square distribution (approximate). More...
 
template<... >
void nonlinearDiffusion (...)
 Perform edge-preserving smoothing at the given scale. More...
 
template<... >
void nonlinearDiffusionExplicit (...)
 Perform edge-preserving smoothing at the given scale using an explicit scheme. More...
 
template<... >
void nonlinearLeastSquares (...)
 Fit a non-linear model to given data by minimizing least squares loss. More...
 
template<... >
bool nonparametricNoiseNormalization (...)
 Noise normalization by means of an estimated non-parametric noise model. More...
 
template<unsigned IntBits, unsigned FracBits>
FixedPoint< IntBits, FracBits > norm (FixedPoint< IntBits, FracBits > const &v)
 norm (same as abs).
 
template<class T >
NormTraits< T >::NormType norm (T const &t)
 The norm of a numerical object. More...
 
template<class R >
FFTWComplex< R >::NormType norm (const FFTWComplex< R > &a)
 norm (= magnitude)
 
template<typename IntType >
Rational< IntType > norm (const Rational< IntType > &r)
 norm (same as abs(r))
 
template<int IntBits, FPOverflowHandling OverflowHandling>
NormTraits< FixedPoint16
< IntBits, OverflowHandling >
>::NormType 
norm (FixedPoint16< IntBits, OverflowHandling > const &v)
 norm (same as abs).
 
template<... >
void normalizedConvolveImage (...)
 Performs a 2-dimensional normalized convolution, i.e. convolution with a mask image. More...
 
template<... >
void normalizedCrossCorrelation (...)
 This function performes a (slow) normalized cross-correlation. More...
 
std::string normalizeString (std::string const &s)
 
bool odd (int t)
 Check if an integer is odd. More...
 
template<class V1 , unsigned int RIDX1, unsigned int GIDX1, unsigned int BIDX1, class V2 , unsigned int RIDX2, unsigned int GIDX2, unsigned int BIDX2>
bool operator!= (RGBValue< V1, RIDX1, GIDX1, BIDX1 > const &l, RGBValue< V2, RIDX2, GIDX2, BIDX2 > const &r)
 component-wise not equal
 
template<unsigned IntBits1, unsigned FracBits1, unsigned IntBits2, unsigned FracBits2>
bool operator!= (FixedPoint< IntBits1, FracBits1 > l, FixedPoint< IntBits2, FracBits2 > r)
 not equal
 
template<class R >
bool operator!= (FFTWComplex< R > const &a, const FFTWComplex< R > &b)
 not equal
 
template<class R >
bool operator!= (FFTWComplex< R > const &a, double b)
 not equal
 
template<class R >
bool operator!= (double a, const FFTWComplex< R > &b)
 not equal
 
template<typename IntType1 , typename IntType2 >
bool operator!= (Rational< IntType1 > const &l, Rational< IntType2 > const &r)
 inequality
 
template<typename IntType1 , typename IntType2 >
bool operator!= (const Rational< IntType1 > &l, IntType2 const &i)
 inequality with right-hand IntType2 argument
 
template<typename IntType1 , typename IntType2 >
bool operator!= (IntType1 const &l, Rational< IntType2 > const &r)
 inequality with left-hand IntType1 argument
 
template<class V1 , int SIZE, class D1 , class D2 , class V2 , class D3 , class D4 >
bool operator!= (TinyVectorBase< V1, SIZE, D1, D2 > const &l, TinyVectorBase< V2, SIZE, D3, D4 > const &r)
 component-wise not equal
 
template<int IntBits1, FPOverflowHandling OverflowHandling, int IntBits2>
bool operator!= (FixedPoint16< IntBits1, OverflowHandling > l, FixedPoint16< IntBits2, OverflowHandling > r)
 not equal
 
template<class V1 , int SIZE, class D1 , class D2 , class V2 , class D3 , class D4 >
PromoteTraits< TinyVector< V1,
SIZE >, TinyVector< V2, SIZE >
>::Promote 
operator% (TinyVectorBase< V1, SIZE, D1, D2 > const &l, TinyVectorBase< V2, SIZE, D3, D4 > const &r)
 component-wise modulo
 
template<typename Type >
Quaternion< Type > operator* (const Quaternion< Type > &t1, const Quaternion< Type > &t2)
 Multiplication.
 
template<typename Type >
Quaternion< Type > operator* (const Quaternion< Type > &t1, double t2)
 Multiplication with a scalar on the right.
 
template<typename Type >
Quaternion< Type > operator* (double t1, const Quaternion< Type > &t2)
 Multiplication with a scalar on the left.
 
template<unsigned IntBits1, unsigned FracBits1, unsigned IntBits2, unsigned FracBits2>
FixedPointTraits< FixedPoint
< IntBits1, FracBits1 >
, FixedPoint< IntBits2,
FracBits2 > >::MultipliesType 
operator* (FixedPoint< IntBits1, FracBits1 > l, FixedPoint< IntBits2, FracBits2 > r)
 multiplication with automatic determination of the appropriate result type.
 
template<class V1 , unsigned int R, unsigned int G, unsigned int B, class V2 >
PromoteTraits< RGBValue< V1, R,
G, B >, RGBValue< V2, R, G, B >
>::Promote 
operator* (RGBValue< V1, R, G, B > const &r1, RGBValue< V2, R, G, B > const &r2)
 component-wise multiplication
 
template<class V , unsigned int R, unsigned int G, unsigned int B>
NumericTraits< RGBValue< V, R,
G, B > >::RealPromote 
operator* (double v, RGBValue< V, R, G, B > const &r)
 component-wise left scalar multiplication
 
template<class V , unsigned int R, unsigned int G, unsigned int B>
NumericTraits< RGBValue< V, R,
G, B > >::RealPromote 
operator* (RGBValue< V, R, G, B > const &r, double v)
 component-wise right scalar multiplication
 
template<typename IntType >
Rational< IntType > operator* (Rational< IntType > l, Rational< IntType > const &r)
 multiplication
 
template<typename IntType >
Rational< IntType > operator* (Rational< IntType > l, typename Rational< IntType >::param_type r)
 multiplication with right-hand IntType argument
 
template<typename IntType >
Rational< IntType > operator* (typename Rational< IntType >::param_type l, Rational< IntType > r)
 multiplication with left-hand IntType argument
 
template<class R >
FFTWComplex< R > operator* (FFTWComplex< R > a, const FFTWComplex< R > &b)
 multiplication
 
template<class R >
FFTWComplex< R > operator* (FFTWComplex< R > a, double b)
 right multiplication with scalar double
 
template<class R >
FFTWComplex< R > operator* (double a, FFTWComplex< R > b)
 left multiplication with scalar double
 
template<int IntBits1, FPOverflowHandling OverflowHandling, int IntBits2>
PromoteTraits< FixedPoint16
< IntBits1, OverflowHandling >
, FixedPoint16< IntBits2,
OverflowHandling > >::Promote 
operator* (FixedPoint16< IntBits1, OverflowHandling > l, FixedPoint16< IntBits2, OverflowHandling > r)
 multiplication with automatic determination of the appropriate result type.
 
template<class V1 , int SIZE, class D1 , class D2 , class V2 , class D3 , class D4 >
PromoteTraits< TinyVector< V1,
SIZE >, TinyVector< V2, SIZE >
>::Promote 
operator* (TinyVectorBase< V1, SIZE, D1, D2 > const &l, TinyVectorBase< V2, SIZE, D3, D4 > const &r)
 component-wise multiplication
 
template<class V , int SIZE, class D1 , class D2 >
NumericTraits< TinyVector< V,
SIZE > >::RealPromote 
operator* (double v, TinyVectorBase< V, SIZE, D1, D2 > const &r)
 component-wise left scalar multiplication
 
template<class V , int SIZE, class D1 , class D2 >
NumericTraits< TinyVector< V,
SIZE > >::RealPromote 
operator* (TinyVectorBase< V, SIZE, D1, D2 > const &l, double v)
 component-wise right scalar multiplication
 
template<class V1 , unsigned int RIDX1, unsigned int GIDX1, unsigned int BIDX1, class V2 , unsigned int RIDX2, unsigned int GIDX2, unsigned int BIDX2>
RGBValue< V1, RIDX1, GIDX1,
BIDX1 > & 
operator*= (RGBValue< V1, RIDX1, GIDX1, BIDX1 > &l, RGBValue< V2, RIDX2, GIDX2, BIDX2 > const &r)
 componentwise multiply-assignment
 
template<class V , unsigned int RIDX, unsigned int GIDX, unsigned int BIDX>
RGBValue< V, RIDX, GIDX, BIDX > & operator*= (RGBValue< V, RIDX, GIDX, BIDX > &l, double r)
 componentwise scalar multiply-assignment
 
template<class R >
FFTWComplex< R > & operator*= (FFTWComplex< R > &a, const FFTWComplex< R > &b)
 multiply-assignment
 
template<class R >
FFTWComplex< R > & operator*= (FFTWComplex< R > &a, double b)
 multiply-assignment with scalar double
 
template<typename Type >
Quaternion< Type > operator+ (const Quaternion< Type > &t1, const Quaternion< Type > &t2)
 Addition.
 
template<typename Type >
Quaternion< Type > operator+ (const Quaternion< Type > &t1, const Type &t2)
 Addition of a scalar on the right.
 
template<typename Type >
Quaternion< Type > operator+ (const Type &t1, const Quaternion< Type > &t2)
 Addition of a scalar on the left.
 
template<unsigned IntBits1, unsigned FracBits1, unsigned IntBits2, unsigned FracBits2>
FixedPointTraits< FixedPoint
< IntBits1, FracBits1 >
, FixedPoint< IntBits2,
FracBits2 > >::PlusType 
operator+ (FixedPoint< IntBits1, FracBits1 > l, FixedPoint< IntBits2, FracBits2 > r)
 addition with automatic determination of the appropriate result type.
 
Diff2D operator+ (Diff2D const &a, Diff2D const &b)
 
Size2D operator+ (Size2D const &a, Diff2D const &b)
 
Point2D operator+ (Point2D const &a, Diff2D const &b)
 
Point2D operator+ (Size2D const &s, Point2D const &p)
 
template<class V1 , unsigned int R, unsigned int G, unsigned int B, class V2 >
PromoteTraits< RGBValue< V1, R,
G, B >, RGBValue< V2, R, G, B >
>::Promote 
operator+ (RGBValue< V1, R, G, B > const &r1, RGBValue< V2, R, G, B > const &r2)
 component-wise addition
 
template<typename IntType >
Rational< IntType > operator+ (const Rational< IntType > &r)
 unary plus
 
template<typename IntType >
Rational< IntType > operator+ (Rational< IntType > l, Rational< IntType > const &r)
 addition
 
template<typename IntType >
Rational< IntType > operator+ (Rational< IntType > l, typename Rational< IntType >::param_type r)
 addition of right-hand IntType argument
 
template<typename IntType >
Rational< IntType > operator+ (typename Rational< IntType >::param_type l, Rational< IntType > r)
 addition of left-hand IntType argument
 
template<class R >
FFTWComplex< R > operator+ (FFTWComplex< R > a, const FFTWComplex< R > &b)
 addition
 
template<class R >
FFTWComplex< R > operator+ (FFTWComplex< R > a, double b)
 right addition with scalar double
 
template<class R >
FFTWComplex< R > operator+ (double a, FFTWComplex< R > b)
 left addition with scalar double
 
template<int IntBits1, FPOverflowHandling OverflowHandling, int IntBits2>
PromoteTraits< FixedPoint16
< IntBits1, OverflowHandling >
, FixedPoint16< IntBits2,
OverflowHandling > >::Promote 
operator+ (FixedPoint16< IntBits1, OverflowHandling > l, FixedPoint16< IntBits2, OverflowHandling > r)
 addition with automatic determination of the appropriate result type.
 
template<class V1 , int SIZE, class D1 , class D2 , class V2 , class D3 , class D4 >
PromoteTraits< TinyVector< V1,
SIZE >, TinyVector< V2, SIZE >
>::Promote 
operator+ (TinyVectorBase< V1, SIZE, D1, D2 > const &l, TinyVectorBase< V2, SIZE, D3, D4 > const &r)
 component-wise addition
 
template<class V , int SIZE, class D1 , class D2 >
NumericTraits< TinyVector< V,
SIZE > >::RealPromote 
operator+ (double v, TinyVectorBase< V, SIZE, D1, D2 > const &r)
 component-wise left scalar addition
 
template<class V , int SIZE, class D1 , class D2 >
NumericTraits< TinyVector< V,
SIZE > >::RealPromote 
operator+ (TinyVectorBase< V, SIZE, D1, D2 > const &l, double v)
 component-wise right scalar addition
 
template<class V1 , unsigned int RIDX1, unsigned int GIDX1, unsigned int BIDX1, class V2 , unsigned int RIDX2, unsigned int GIDX2, unsigned int BIDX2>
RGBValue< V1, RIDX1, GIDX1,
BIDX1 > & 
operator+= (RGBValue< V1, RIDX1, GIDX1, BIDX1 > &l, RGBValue< V2, RIDX2, GIDX2, BIDX2 > const &r)
 componentwise add-assignment
 
template<class R >
FFTWComplex< R > & operator+= (FFTWComplex< R > &a, const FFTWComplex< R > &b)
 add-assignment
 
template<class R >
FFTWComplex< R > & operator+= (FFTWComplex< R > &a, double b)
 add-assignment with scalar double
 
template<typename Type >
Quaternion< Type > operator- (const Quaternion< Type > &t1, const Quaternion< Type > &t2)
 Subtraction.
 
template<typename Type >
Quaternion< Type > operator- (const Quaternion< Type > &t1, const Type &t2)
 Subtraction of a scalar on the right.
 
template<typename Type >
Quaternion< Type > operator- (const Type &t1, const Quaternion< Type > &t2)
 Subtraction of a scalar on the left.
 
template<unsigned IntBits1, unsigned FracBits1, unsigned IntBits2, unsigned FracBits2>
FixedPointTraits< FixedPoint
< IntBits1, FracBits1 >
, FixedPoint< IntBits2,
FracBits2 > >::MinusType 
operator- (FixedPoint< IntBits1, FracBits1 > l, FixedPoint< IntBits2, FracBits2 > r)
 subtraction with automatic determination of the appropriate result type.
 
Diff2D operator- (Diff2D const &a, Diff2D const &b)
 
Size2D operator- (Size2D const &s, Diff2D const &offset)
 
Point2D operator- (Point2D const &s, Diff2D const &offset)
 
Size2D operator- (Point2D const &s, Point2D const &p)
 
template<class V1 , unsigned int R, unsigned int G, unsigned int B, class V2 >
PromoteTraits< RGBValue< V1, R,
G, B >, RGBValue< V2, R, G, B >
>::Promote 
operator- (RGBValue< V1, R, G, B > const &r1, RGBValue< V2, R, G, B > const &r2)
 component-wise subtraction
 
template<typename IntType >
Rational< IntType > operator- (const Rational< IntType > &r)
 unary minus (negation)
 
template<typename IntType >
Rational< IntType > operator- (Rational< IntType > l, Rational< IntType > const &r)
 subtraction
 
template<typename IntType >
Rational< IntType > operator- (Rational< IntType > l, typename Rational< IntType >::param_type r)
 subtraction of right-hand IntType argument
 
template<typename IntType >
Rational< IntType > operator- (typename Rational< IntType >::param_type l, Rational< IntType > const &r)
 subtraction from left-hand IntType argument
 
template<class R >
FFTWComplex< R > operator- (FFTWComplex< R > a, const FFTWComplex< R > &b)
 subtraction
 
template<class R >
FFTWComplex< R > operator- (FFTWComplex< R > a, double b)
 right subtraction with scalar double
 
template<class R >
FFTWComplex< R > operator- (double a, FFTWComplex< R > const &b)
 left subtraction with scalar double
 
template<int IntBits1, FPOverflowHandling OverflowHandling, int IntBits2>
PromoteTraits< FixedPoint16
< IntBits1, OverflowHandling >
, FixedPoint16< IntBits2,
OverflowHandling > >::Promote 
operator- (FixedPoint16< IntBits1, OverflowHandling > l, FixedPoint16< IntBits2, OverflowHandling > r)
 subtraction with automatic determination of the appropriate result type.
 
template<class V1 , int SIZE, class D1 , class D2 , class V2 , class D3 , class D4 >
PromoteTraits< TinyVector< V1,
SIZE >, TinyVector< V2, SIZE >
>::Promote 
operator- (TinyVectorBase< V1, SIZE, D1, D2 > const &l, TinyVectorBase< V2, SIZE, D3, D4 > const &r)
 component-wise subtraction
 
template<class V , int SIZE, class D1 , class D2 >
NumericTraits< TinyVector< V,
SIZE > >::RealPromote 
operator- (double v, TinyVectorBase< V, SIZE, D1, D2 > const &r)
 component-wise left scalar subtraction
 
template<class V , int SIZE, class D1 , class D2 >
NumericTraits< TinyVector< V,
SIZE > >::RealPromote 
operator- (TinyVectorBase< V, SIZE, D1, D2 > const &l, double v)
 component-wise right scalar subtraction
 
template<class V , int SIZE, class D1 , class D2 >
TinyVector< V, SIZE > operator- (TinyVectorBase< V, SIZE, D1, D2 > const &v)
 
template<class V1 , unsigned int RIDX1, unsigned int GIDX1, unsigned int BIDX1, class V2 , unsigned int RIDX2, unsigned int GIDX2, unsigned int BIDX2>
RGBValue< V1, RIDX1, GIDX1,
BIDX1 > & 
operator-= (RGBValue< V1, RIDX1, GIDX1, BIDX1 > &l, RGBValue< V2, RIDX2, GIDX2, BIDX2 > const &r)
 componentwise subtract-assignment
 
template<class R >
FFTWComplex< R > & operator-= (FFTWComplex< R > &a, const FFTWComplex< R > &b)
 subtract-assignment
 
template<class R >
FFTWComplex< R > & operator-= (FFTWComplex< R > &a, double b)
 subtract-assignment with scalar double
 
template<typename Type >
Quaternion< Type > operator/ (const Quaternion< Type > &t1, const Quaternion< Type > &t2)
 Division.
 
template<typename Type >
Quaternion< Type > operator/ (const Quaternion< Type > &t1, double t2)
 Division by a scalar.
 
template<typename Type >
Quaternion< Type > operator/ (double t1, const Quaternion< Type > &t2)
 Division of a scalar by a Quaternion.
 
template<class V1 , unsigned int R, unsigned int G, unsigned int B, class V2 >
PromoteTraits< RGBValue< V1, R,
G, B >, RGBValue< V2, R, G, B >
>::Promote 
operator/ (RGBValue< V1, R, G, B > const &r1, RGBValue< V2, R, G, B > const &r2)
 component-wise division
 
template<typename IntType >
Rational< IntType > operator/ (Rational< IntType > l, Rational< IntType > const &r)
 division
 
template<class V , unsigned int R, unsigned int G, unsigned int B>
NumericTraits< RGBValue< V, R,
G, B > >::RealPromote 
operator/ (RGBValue< V, R, G, B > const &r, double v)
 component-wise scalar division
 
template<typename IntType >
Rational< IntType > operator/ (Rational< IntType > l, typename Rational< IntType >::param_type r)
 division by right-hand IntType argument
 
template<typename IntType >
Rational< IntType > operator/ (typename Rational< IntType >::param_type l, Rational< IntType > const &r)
 division of left-hand IntType argument
 
template<class R >
FFTWComplex< R > operator/ (FFTWComplex< R > a, const FFTWComplex< R > &b)
 division
 
template<class R >
FFTWComplex< R > operator/ (FFTWComplex< R > a, double b)
 right division with scalar double
 
template<int IntBits1, FPOverflowHandling OverflowHandling, int IntBits2>
PromoteTraits< FixedPoint16
< IntBits1, OverflowHandling >
, FixedPoint16< IntBits2,
OverflowHandling > >::Promote 
operator/ (FixedPoint16< IntBits1, OverflowHandling > l, FixedPoint16< IntBits2, OverflowHandling > r)
 division with automatic determination of the appropriate result type.
 
template<class V1 , int SIZE, class D1 , class D2 , class V2 , class D3 , class D4 >
PromoteTraits< TinyVector< V1,
SIZE >, TinyVector< V2, SIZE >
>::Promote 
operator/ (TinyVectorBase< V1, SIZE, D1, D2 > const &l, TinyVectorBase< V2, SIZE, D3, D4 > const &r)
 component-wise division
 
template<class V , int SIZE, class D1 , class D2 >
NumericTraits< TinyVector< V,
SIZE > >::RealPromote 
operator/ (double v, TinyVectorBase< V, SIZE, D1, D2 > const &r)
 component-wise left scalar division
 
template<class V , int SIZE, class D1 , class D2 >
NumericTraits< TinyVector< V,
SIZE > >::RealPromote 
operator/ (TinyVectorBase< V, SIZE, D1, D2 > const &l, double v)
 component-wise right scalar division
 
template<class V1 , unsigned int RIDX1, unsigned int GIDX1, unsigned int BIDX1, class V2 , unsigned int RIDX2, unsigned int GIDX2, unsigned int BIDX2>
RGBValue< V1, RIDX1, GIDX1,
BIDX1 > & 
operator/= (RGBValue< V1, RIDX1, GIDX1, BIDX1 > &l, RGBValue< V2, RIDX2, GIDX2, BIDX2 > const &r)
 componentwise divide-assignment
 
template<class V , unsigned int RIDX, unsigned int GIDX, unsigned int BIDX>
RGBValue< V, RIDX, GIDX, BIDX > & operator/= (RGBValue< V, RIDX, GIDX, BIDX > &l, double r)
 componentwise scalar divide-assignment
 
template<class R >
FFTWComplex< R > & operator/= (FFTWComplex< R > &a, const FFTWComplex< R > &b)
 divide-assignment
 
template<class R >
FFTWComplex< R > & operator/= (FFTWComplex< R > &a, double b)
 divide-assignment with scalar double
 
template<unsigned IntBits1, unsigned FracBits1, unsigned IntBits2, unsigned FracBits2>
bool operator< (FixedPoint< IntBits1, FracBits1 > l, FixedPoint< IntBits2, FracBits2 > r)
 less than
 
template<typename IntType1 , typename IntType2 >
bool operator< (const Rational< IntType1 > &l, const Rational< IntType2 > &r)
 less-than
 
template<typename IntType1 , typename IntType2 >
bool operator< (const Rational< IntType1 > &l, IntType2 const &i)
 less-than with right-hand IntType2 argument
 
template<typename IntType1 , typename IntType2 >
bool operator< (IntType1 const &l, Rational< IntType2 > const &r)
 less-than with left-hand IntType1 argument
 
template<class V1 , int SIZE, class D1 , class D2 , class V2 , class D3 , class D4 >
bool operator< (TinyVectorBase< V1, SIZE, D1, D2 > const &l, TinyVectorBase< V2, SIZE, D3, D4 > const &r)
 lexicographical comparison
 
template<int IntBits1, FPOverflowHandling OverflowHandling, int IntBits2>
bool operator< (FixedPoint16< IntBits1, OverflowHandling > l, FixedPoint16< IntBits2, OverflowHandling > r)
 less than
 
std::ostream & operator<< (std::ostream &os, const SIFImportInfo &info)
 
template<class V1 , int SIZE, class DATA , class DERIVED >
std::ostream & operator<< (std::ostream &out, TinyVectorBase< V1, SIZE, DATA, DERIVED > const &l)
 stream output
 
template<unsigned IntBits1, unsigned FracBits1, unsigned IntBits2, unsigned FracBits2>
bool operator<= (FixedPoint< IntBits1, FracBits1 > l, FixedPoint< IntBits2, FracBits2 > r)
 less or equal
 
template<typename IntType1 , typename IntType2 >
bool operator<= (Rational< IntType1 > const &l, Rational< IntType2 > const &r)
 less-equal
 
template<typename IntType1 , typename IntType2 >
bool operator<= (Rational< IntType1 > const &l, IntType2 const &r)
 less-equal with right-hand IntType2 argument
 
template<typename IntType1 , typename IntType2 >
bool operator<= (IntType1 const &l, Rational< IntType2 > const &r)
 less-equal with left-hand IntType1 argument
 
template<int IntBits1, FPOverflowHandling OverflowHandling, int IntBits2>
bool operator<= (FixedPoint16< IntBits1, OverflowHandling > l, FixedPoint16< IntBits2, OverflowHandling > r)
 less or equal
 
template<class V1 , unsigned int RIDX1, unsigned int GIDX1, unsigned int BIDX1, class V2 , unsigned int RIDX2, unsigned int GIDX2, unsigned int BIDX2>
bool operator== (RGBValue< V1, RIDX1, GIDX1, BIDX1 > const &l, RGBValue< V2, RIDX2, GIDX2, BIDX2 > const &r)
 component-wise equal
 
template<unsigned IntBits1, unsigned FracBits1, unsigned IntBits2, unsigned FracBits2>
bool operator== (FixedPoint< IntBits1, FracBits1 > l, FixedPoint< IntBits2, FracBits2 > r)
 equal
 
template<class R >
bool operator== (FFTWComplex< R > const &a, const FFTWComplex< R > &b)
 equal
 
template<typename IntType1 , typename IntType2 >
bool operator== (const Rational< IntType1 > &l, const Rational< IntType2 > &r)
 equality
 
template<typename IntType1 , typename IntType2 >
bool operator== (const Rational< IntType1 > &l, IntType2 const &i)
 equality with right-hand IntType2 argument
 
template<typename IntType1 , typename IntType2 >
bool operator== (IntType1 const &l, Rational< IntType2 > const &r)
 equality with left-hand IntType1 argument
 
template<class V1 , int SIZE, class D1 , class D2 , class V2 , class D3 , class D4 >
bool operator== (TinyVectorBase< V1, SIZE, D1, D2 > const &l, TinyVectorBase< V2, SIZE, D3, D4 > const &r)
 component-wise equal
 
template<int IntBits1, FPOverflowHandling OverflowHandling, int IntBits2>
bool operator== (FixedPoint16< IntBits1, OverflowHandling > l, FixedPoint16< IntBits2, OverflowHandling > r)
 equal
 
template<unsigned IntBits1, unsigned FracBits1, unsigned IntBits2, unsigned FracBits2>
bool operator> (FixedPoint< IntBits1, FracBits1 > l, FixedPoint< IntBits2, FracBits2 > r)
 greater
 
template<typename IntType1 , typename IntType2 >
bool operator> (Rational< IntType1 > const &l, Rational< IntType2 > const &r)
 greater-than
 
template<typename IntType1 , typename IntType2 >
bool operator> (const Rational< IntType1 > &l, IntType2 const &i)
 greater-than with right-hand IntType2 argument
 
template<typename IntType1 , typename IntType2 >
bool operator> (IntType1 const &l, Rational< IntType2 > const &r)
 greater-than with left-hand IntType1 argument
 
template<int IntBits1, FPOverflowHandling OverflowHandling, int IntBits2>
bool operator> (FixedPoint16< IntBits1, OverflowHandling > l, FixedPoint16< IntBits2, OverflowHandling > r)
 greater
 
template<unsigned IntBits1, unsigned FracBits1, unsigned IntBits2, unsigned FracBits2>
bool operator>= (FixedPoint< IntBits1, FracBits1 > l, FixedPoint< IntBits2, FracBits2 > r)
 greater or equal
 
template<typename IntType1 , typename IntType2 >
bool operator>= (Rational< IntType1 > const &l, Rational< IntType2 > const &r)
 greater-equal
 
template<typename IntType1 , typename IntType2 >
bool operator>= (Rational< IntType1 > const &l, IntType2 const &r)
 greater-equal with right-hand IntType2 argument
 
template<typename IntType1 , typename IntType2 >
bool operator>= (IntType1 const &l, Rational< IntType2 > const &r)
 greater-equal with left-hand IntType1 argument
 
template<int IntBits1, FPOverflowHandling OverflowHandling, int IntBits2>
bool operator>= (FixedPoint16< IntBits1, OverflowHandling > l, FixedPoint16< IntBits2, OverflowHandling > r)
 greater or equal
 
template<... >
void parallel_foreach (...)
 Apply a functor to all items in a range in parallel. More...
 
template<class NODE , class PREDECESSORS >
size_t pathLength (const NODE source, const NODE target, const PREDECESSORS &predecessors)
 get the length in node units of a path
 
template<... >
void pLSA (...)
 Decompose a matrix according to the pLSA algorithm. More...
 
TinyVector< float, 3 > polar2Lab (double color, double brightness, double saturation)
 Init L*a*b* color triple from polar representation. More...
 
TinyVector< float, 3 > polar2Luv (double color, double brightness, double saturation)
 Init L*u*v* color triple from polar representation. More...
 
TinyVector< float, 3 > polar2YPrimeCbCr (double color, double brightness, double saturation)
 Init Y'CbCr color triple from polar representation. More...
 
TinyVector< float, 3 > polar2YPrimeIQ (double color, double brightness, double saturation)
 Init Y'IQ color triple from polar representation. More...
 
TinyVector< float, 3 > polar2YPrimePbPr (double color, double brightness, double saturation)
 Init Y'PbPr color triple from polar representation. More...
 
TinyVector< float, 3 > polar2YPrimeUV (double color, double brightness, double saturation)
 Init Y'UV color triple from polar representation. More...
 
template<int PolynomOrder, class SrcPointIterator , class DestPointIterator >
linalg::TemporaryMatrix< double > polynomialMatrix2DFromCorrespondingPoints (SrcPointIterator s, SrcPointIterator s_end, DestPointIterator d)
 Create polynomial matrix of a certain degree that maps corresponding points onto each other. More...
 
template<class POLYNOMIAL , class VECTOR >
bool polynomialRealRoots (POLYNOMIAL const &p, VECTOR &roots, bool polishRoots)
 
template<class POLYNOMIAL , class VECTOR >
bool polynomialRoots (POLYNOMIAL const &poriginal, VECTOR &roots, bool polishRoots)
 
template<... >
void polynomialWarpImage (...)
 Warp an image according to an polynomial transformation. More...
 
std::vector< double > polynomialWarpWeights (double x, double y, unsigned int polynom_order)
 
template<typename IntType >
Rational< IntType > pow (const Rational< IntType > &r, int n)
 
template<class V , int SIZE, class D1 , class D2 , class E >
TinyVector< V, SIZE > pow (TinyVectorBase< V, SIZE, D1, D2 > const &v, E exponent)
 
template<unsigned n, class V >
power (const V &x)
 Exponentiation to a positive integer power by squaring. More...
 
template<class T , class C1 , class C2 , class C3 >
void principalComponents (MultiArrayView< 2, T, C1 > const &features, MultiArrayView< 2, T, C2 > fz, MultiArrayView< 2, T, C3 > zv)
 Decompose a matrix according to the PCA algorithm. More...
 
template<class V , int SIZE, class D1 , class D2 >
NumericTraits< V >::Promote prod (TinyVectorBase< V, SIZE, D1, D2 > const &l)
 product of the vector's elements
 
template<class BASE_GRAPH , class BASE_GRAPH_LABELS , class RAG_FEATURES , class BASE_GRAPH_FEATURES >
void projectBack (const AdjacencyListGraph &rag, const BASE_GRAPH &bg, const Int64 ignoreLabel, const BASE_GRAPH_LABELS bgLabels, const RAG_FEATURES &ragFeatures, BASE_GRAPH_FEATURES &bgFeatures)
 
template<class RAG , class BASE_GRAPH , class BASE_GRAPH_RAG_LABELS , class BASE_GRAPH_GT , class RAG_GT , class RAG_GT_QT >
void projectGroundTruth (const RAG &rag, const BASE_GRAPH &baseGraph, const BASE_GRAPH_RAG_LABELS &baseGraphRagLabels, const BASE_GRAPH_GT &baseGraphGt, RAG_GT &ragGt, RAG_GT_QT &)
 
template<class SrcPointIterator , class DestPointIterator >
linalg::TemporaryMatrix< double > projectiveMatrix2DFromCorrespondingPoints (SrcPointIterator s, SrcPointIterator send, DestPointIterator d)
 Create homogeneous matrix that maps corresponding points onto each other. More...
 
template<... >
void projectiveWarpImage (...)
 Warp an image according to an projective transformation. More...
 
template<... >
void pyramidExpandBurtFilter (...)
 Two-fold up-sampling for image pyramid reconstruction. More...
 
template<class Image , class Alloc >
void pyramidExpandBurtLaplacian (ImagePyramid< Image, Alloc > &pyramid, int fromLevel, int toLevel, double centerValue=0.4)
 Reconstruct a Laplacian pyramid. More...
 
template<... >
void pyramidReduceBurtFilter (...)
 Two-fold down-sampling for image pyramid construction. More...
 
template<class Image , class Alloc >
void pyramidReduceBurtLaplacian (ImagePyramid< Image, Alloc > &pyramid, int fromLevel, int toLevel, double centerValue=0.4)
 Create a Laplacian pyramid. More...
 
template<... >
bool quadraticNoiseNormalization (...)
 Noise normalization by means of an estimated or given quadratic noise model. More...
 
template<... >
unsigned int quadraticProgramming (...)
 
double radialGaborSigma (double centerFrequency)
 Calculate sensible radial sigma for given parameters. More...
 
template<... >
void radialSymmetryTransform (...)
 Find centers of radial symmetry in an image. More...
 
template<class RAGGRAPH , class GRAPH , class RAGEDGES , unsigned int N, class T >
MultiArray< 2, MultiArrayIndexragFindEdges (const RAGGRAPH &rag, const GRAPH &graph, const RAGEDGES &affiliatedEdges, MultiArrayView< N, T > labelsArray, const typename RAGGRAPH::Node &node)
 Find indices of points on the edges. More...
 
RandomMT19937randomMT19937 ()
 
RandomTT800randomTT800 ()
 
template<typename T , typename IntType >
rational_cast (const Rational< IntType > &src)
 
template<class RadialBasisFunctor , class SrcPointIterator , class DestPointIterator >
linalg::TemporaryMatrix< double > rbfMatrix2DFromCorrespondingPoints (SrcPointIterator s, SrcPointIterator s_end, DestPointIterator d, RadialBasisFunctor const &rbf)
 Create a matrix that maps corresponding points onto each other using a given RBF. More...
 
template<int ORDER, class T , class DestIterator , class DestAccessor , class DestPointIterator , class C , class RadialBasisFunctor >
void rbfWarpImage (SplineImageView< ORDER, T > const &src, DestIterator dul, DestIterator dlr, DestAccessor dest, DestPointIterator d, DestPointIterator d_end, MultiArrayView< 2, double, C > const &W, RadialBasisFunctor rbf)
 Warp an image according to an radial basis function based transformation. More...
 
template<... >
void readHDF5 (...)
 Read the data specified by the given vigra::HDF5ImportInfo object and write the into the given 'array'. More...
 
void readSIF (const SIFImportInfo &info, MultiArrayView< 3, float > array)
 Read the image data specified by the given vigra::SIFImportInfo object and write them into the given 'array'. More...
 
void readSIFBlock (const SIFImportInfo &info, Shape3 offset, Shape3 shape, MultiArrayView< 3, float > array)
 Read parts of the image data from an Andor SIF file specified with an SIFImportInfo object and write them into the MultiArray array. More...
 
template<class R >
real (const FFTWComplex< R > &a)
 real part
 
template<... >
void recursiveFilterLine (...)
 Performs a 1-dimensional recursive convolution of the source signal. More...
 
template<... >
void recursiveFilterX (...)
 Performs 1 dimensional recursive filtering (1st and 2nd order) in x direction. More...
 
template<... >
void recursiveFilterY (...)
 Performs 1 dimensional recursive filtering (1st and 2nd order) in y direction. More...
 
template<... >
void recursiveFirstDerivativeLine (...)
 Performs a 1 dimensional recursive convolution of the source signal. More...
 
template<... >
void recursiveFirstDerivativeX (...)
 Recursively calculates the 1 dimensional first derivative in x direction. More...
 
template<... >
void recursiveFirstDerivativeY (...)
 Recursively calculates the 1 dimensional first derivative in y direction. More...
 
template<... >
void recursiveGaussianFilterLine (...)
 Compute a 1-dimensional recursive approximation of Gaussian smoothing. More...
 
template<... >
void recursiveGaussianFilterX (...)
 Compute 1 dimensional recursive approximation of Gaussian smoothing in y direction. More...
 
template<... >
void recursiveGaussianFilterY (...)
 Compute 1 dimensional recursive approximation of Gaussian smoothing in y direction. More...
 
template<class GRAPH , class NODE_FEATURES_IN , class EDGE_INDICATOR , class NODE_FEATURES_OUT >
void recursiveGraphSmoothing (const GRAPH &g, const NODE_FEATURES_IN &nodeFeaturesIn, const EDGE_INDICATOR &edgeIndicator, const float lambda, const float edgeThreshold, const float scale, size_t iterations, NODE_FEATURES_OUT &nodeFeaturesBuffer, NODE_FEATURES_OUT &nodeFeaturesOut)
 smooth node features of a graph More...
 
template<... >
void recursiveSecondDerivativeLine (...)
 Performs a 1 dimensional recursive convolution of the source signal. More...
 
template<... >
void recursiveSecondDerivativeX (...)
 Recursively calculates the 1 dimensional second derivative in x direction. More...
 
template<... >
void recursiveSecondDerivativeY (...)
 Recursively calculates the 1 dimensional second derivative in y direction. More...
 
template<... >
void recursiveSmoothLine (...)
 Convolves the image with a 1-dimensional exponential filter. More...
 
template<... >
void recursiveSmoothX (...)
 Performs 1 dimensional recursive smoothing in x direction. More...
 
template<... >
void recursiveSmoothY (...)
 Performs 1 dimensional recursive smoothing in y direction. More...
 
template<... >
void reflectImage (...)
 Reflect image horizontally or vertically. More...
 
template<... >
void regionImageToCrackEdgeImage (...)
 Transform a labeled image into a crack edge (interpixel edge) image. More...
 
template<... >
void regionImageToEdgeImage (...)
 Transform a labeled image into an edge image. More...
 
template<... >
void removeShortEdges (...)
 Remove short edges from an edge image. More...
 
template<... >
void resampleImage (...)
 Resample image by a given factor. More...
 
template<... >
void resamplingConvolveImage (...)
 Apply two separable resampling filters successively, the first in x-direction, the second in y-direction. More...
 
template<... >
void resamplingConvolveLine (...)
 Performs a 1-dimensional resampling convolution of the source signal using the given set of kernels. More...
 
template<... >
void resamplingConvolveX (...)
 Apply a resampling filter in the x-direction. More...
 
template<... >
void resamplingConvolveY (...)
 Apply a resampling filter in the y-direction. More...
 
template<... >
void resizeImageCatmullRomInterpolation (...)
 Resize image using the Catmull/Rom interpolation function. More...
 
template<... >
void resizeImageCoscotInterpolation (...)
 Resize image using the Coscot interpolation function. More...
 
template<... >
void resizeImageLinearInterpolation (...)
 Resize image using linear interpolation. More...
 
template<... >
void resizeImageNoInterpolation (...)
 Resize image by repeating the nearest pixel values. More...
 
template<... >
void resizeImageSplineInterpolation (...)
 Resize image using B-spline interpolation. More...
 
template<... >
void resizeMultiArraySplineInterpolation (...)
 Resize MultiArray using B-spline interpolation. More...
 
template<class V , int SIZE>
TinyVector< V, SIZE > reverse (TinyVector< V, SIZE > const &t)
 reversed copy
 
detail::RF_DEFAULT & rf_default ()
 factory function to return a RF_DEFAULT tag More...
 
template<class T , class Tag >
void rf_export_HDF5 (const RandomForest< T, Tag > &rf, HDF5File &h5context, const std::string &pathname="")
 Save a random forest to an HDF5File object into a specified HDF5 group. More...
 
template<class T , class Tag >
void rf_export_HDF5 (const RandomForest< T, Tag > &rf, const std::string &filename, const std::string &pathname="")
 Save a random forest to a named HDF5 file into a specified HDF5 group. More...
 
template<class T , class Tag >
void rf_export_HDF5 (const RandomForest< T, Tag > &rf, hid_t outf_id, const std::string &pathname="")
 Save a random forest to an HDF5 file specified by its id. More...
 
template<class T , class Tag >
bool rf_import_HDF5 (RandomForest< T, Tag > &rf, HDF5File &h5context, const std::string &pathname="")
 Read a random forest from an HDF5File object's specified group. More...
 
template<class T , class Tag >
bool rf_import_HDF5 (RandomForest< T, Tag > &rf, const std::string &filename, const std::string &pathname="")
 Read a random forest from a named HDF5 file's specified group. More...
 
template<class T , class Tag >
bool rf_import_HDF5 (RandomForest< T, Tag > &rf, hid_t inf_id, const std::string &pathname="")
 Read a random forest from an HDF5 file specified by its id. More...
 
template<... >
void rieszTransformOfLOG (...)
 Calculate Riesz transforms of the Laplacian of Gaussian. More...
 
template<... >
void rohrCornerDetector (...)
 Find corners in an image (3). More...
 
template<... >
void rotateImage (...)
 Rotate an image by a multiple of 90 degrees or by an arbitrary angle. More...
 
linalg::TemporaryMatrix< double > rotationMatrix2DDegrees (double angle)
 Create homogeneous matrix representing a 2D rotation about the coordinate origin. More...
 
linalg::TemporaryMatrix< double > rotationMatrix2DDegrees (double angle, TinyVector< double, 2 > const &center)
 Create homogeneous matrix representing a 2D rotation about the given point. More...
 
linalg::TemporaryMatrix< double > rotationMatrix2DRadians (double angle)
 Create homogeneous matrix representing a 2D rotation about the coordinate origin. More...
 
linalg::TemporaryMatrix< double > rotationMatrix2DRadians (double angle, TinyVector< double, 2 > const &center)
 Create homogeneous matrix representing a 2D rotation about the given point. More...
 
REAL round (REAL v)
 The rounding function. More...
 
template<unsigned IntBits, unsigned FracBits>
int round (FixedPoint< IntBits, FracBits > v)
 rounding to the nearest integer.
 
template<int IntBits, FPOverflowHandling OverflowHandling>
Int32 round (FixedPoint16< IntBits, OverflowHandling > v)
 rounding to the nearest integer.
 
template<class V , int SIZE, class D1 , class D2 >
TinyVector< V, SIZE > round (TinyVectorBase< V, SIZE, D1, D2 > const &v)
 
long long roundi (double t)
 Round and cast to integer. More...
 
template<int IntBits, FPOverflowHandling OverflowHandling>
Int32 roundi (FixedPoint16< IntBits, OverflowHandling > v)
 rounding to the nearest integer.
 
template<class V , int SIZE, class D1 , class D2 >
TinyVector< std::ptrdiff_t, SIZE > roundi (TinyVectorBase< V, SIZE, D1, D2 > const &v)
 
linalg::TemporaryMatrix< double > scalingMatrix2D (double scalingFactor)
 Create homogeneous matrix representing a 2D uniform scaling about the coordinate origin. More...
 
linalg::TemporaryMatrix< double > scalingMatrix2D (double sx, double sy)
 Create homogeneous matrix representing a 2D non-uniform scaling about the coordinate origin. More...
 
template<... >
void secondOrderTotalVariationFilter (...)
 Performs Anisotropic Total Variation Regularization. More...
 
template<... >
void seededRegionGrowing (...)
 Region Segmentation by means of Seeded Region Growing. More...
 
template<... >
void seededRegionGrowing3D (...)
 Three-dimensional Region Segmentation by means of Seeded Region Growing. More...
 
template<... >
void separableConvolveMultiArray (...)
 Separated convolution on multi-dimensional arrays. More...
 
template<... >
void separableConvolveX (...)
 Performs a 1 dimensional convolution in x direction. More...
 
template<... >
void separableConvolveY (...)
 Performs a 1 dimensional convolution in y direction. More...
 
template<... >
void separableMultiDistance (...)
 Euclidean distance on multi-dimensional arrays. More...
 
template<... >
void separableMultiDistSquared (...)
 Euclidean distance squared on multi-dimensional arrays. More...
 
template<... >
void separableVectorDistance (...)
 Compute the vector distance transform of a N-dimensional binary array. More...
 
linalg::TemporaryMatrix< double > shearMatrix2D (double s01, double s10)
 Create homogeneous matrix representing a 2D shearing. More...
 
template<class GRAPH , class WEIGHTS , class PREDECESSORS , class DISTANCE , class HEURSTIC >
void shortestPathAStar (const GRAPH &graph, const typename GRAPH::Node &source, const typename GRAPH::Node &target, const WEIGHTS &weights, PREDECESSORS &predecessors, DISTANCE &distance, const HEURSTIC &heuristic)
 Astar Shortest path search.
 
template<class T >
sign (T t)
 The sign function. More...
 
template<class T1 , class T2 >
T1 sign (T1 t1, T2 t2)
 The binary sign function. More...
 
template<class T >
int signi (T t)
 The integer sign function. More...
 
template<... >
void simpleSharpening (...)
 Perform simple sharpening function. More...
 
template<class REAL >
REAL sin_pi (REAL x)
 sin(pi*x). More...
 
template<... >
void skeletonizeImage (...)
 Skeletonization of all regions in a labeled 2D image. More...
 
template<... >
void slantedEdgeMTF (...)
 Determine the magnitude transfer function of the camera. More...
 
template<... >
unsigned int slicSuperpixels (...)
 Compute SLIC superpixels in arbitrary dimensions. More...
 
template<class T >
NumericTraits< T >::Promote sq (T t)
 The square function. More...
 
template<unsigned IntBits, unsigned FracBits>
SquareRootTraits< FixedPoint
< IntBits, FracBits >
>::SquareRootResult 
sqrt (FixedPoint< IntBits, FracBits > v)
 square root.
 
template<int IntBits, FPOverflowHandling OverflowHandling>
SquareRootTraits< FixedPoint16
< IntBits, OverflowHandling >
>::SquareRootResult 
sqrt (FixedPoint16< IntBits, OverflowHandling > v)
 square root.
 
template<class V , int SIZE, class D1 , class D2 >
TinyVector< V, SIZE > sqrt (TinyVectorBase< V, SIZE, D1, D2 > const &v)
 
Int32 sqrti (Int32 v)
 Signed integer square root. More...
 
template<class V1 , int SIZE, class D1 , class D2 >
TinyVectorBase< V1, SIZE, D1,
D2 >::SquaredNormType 
squaredNorm (TinyVectorBase< V1, SIZE, D1, D2 > const &t)
 squared norm
 
template<typename Type >
Quaternion< Type >::SquaredNormType squaredNorm (Quaternion< Type > const &q)
 squared norm
 
template<unsigned IntBits, unsigned FracBits>
FixedPointTraits< FixedPoint
< IntBits, FracBits >
, FixedPoint< IntBits,
FracBits > >::MultipliesType 
squaredNorm (FixedPoint< IntBits, FracBits > v)
 squared norm (same as v*v).
 
NormTraits< T >::SquaredNormType squaredNorm (T const &t)
 The squared norm of a numerical object. More...
 
template<class R >
FFTWComplex< R >::SquaredNormType squaredNorm (const FFTWComplex< R > &a)
 squared norm (= squared magnitude)
 
template<typename IntType >
NormTraits< Rational< IntType >
>::SquaredNormType 
squaredNorm (const Rational< IntType > &r)
 squared norm
 
template<int IntBits, FPOverflowHandling OverflowHandling>
NormTraits< FixedPoint16
< IntBits, OverflowHandling >
>::SquaredNormType 
squaredNorm (FixedPoint16< IntBits, OverflowHandling > v)
 squared norm (same as v*v).
 
template<class V , int SIZE>
TinyVector< V, SIZE >
::SquaredNormType 
squaredNorm (TinyVector< V, SIZE > const &t)
 squared norm
 
template<... >
void structureTensor (...)
 Calculate the Structure Tensor for each pixel of and image, using Gaussian (derivative) filters. More...
 
template<... >
void structureTensorMultiArray (...)
 Calculate the structure tensor of a multi-dimensional arrays. More...
 
template<unsigned IntBits1, unsigned FracBits1, unsigned IntBits2, unsigned FracBits2, unsigned IntBits3, unsigned FracBits3>
void sub (FixedPoint< IntBits1, FracBits1 > l, FixedPoint< IntBits2, FracBits2 > r, FixedPoint< IntBits3, FracBits3 > &result)
 subtraction with enforced result type.
 
template<int IntBits1, FPOverflowHandling OverflowHandling, int IntBits2, int IntBits3>
FixedPoint16< IntBits3,
OverflowHandling > & 
sub (FixedPoint16< IntBits1, OverflowHandling > l, FixedPoint16< IntBits2, OverflowHandling > r, FixedPoint16< IntBits3, OverflowHandling > &result)
 subtraction with enforced result type.
 
template<class V , int SIZE, class D1 , class D2 >
NumericTraits< V >::Promote sum (TinyVectorBase< V, SIZE, D1, D2 > const &l)
 sum of the vector's elements
 
template<class T >
void symmetric2x2Eigenvalues (T a00, T a01, T a11, T *r0, T *r1)
 Compute the eigenvalues of a 2x2 real symmetric matrix. More...
 
template<class T >
void symmetric3x3Eigenvalues (T a00, T a01, T a02, T a11, T a12, T a22, T *r0, T *r1, T *r2)
 Compute the eigenvalues of a 3x3 real symmetric matrix. More...
 
template<... >
void symmetricGradientMultiArray (...)
 Calculate gradient of a multi-dimensional arrays using symmetric difference filters. More...
 
template<... >
void tensorDeterminantMultiArray (...)
 Calculate the tensor determinant for every element of a ND tensor array. More...
 
template<... >
void tensorEigenRepresentation (...)
 Calculate eigen representation of a symmetric 2x2 tensor. More...
 
template<... >
void tensorEigenvaluesMultiArray (...)
 Calculate the tensor eigenvalues for every element of a N-D tensor array. More...
 
template<... >
void tensorToEdgeCorner (...)
 Decompose a symmetric 2x2 tensor into its edge and corner parts. More...
 
template<... >
void tensorTrace (...)
 Calculate the trace of a 2x2 tensor. More...
 
template<... >
void tensorTraceMultiArray (...)
 Calculate the tensor trace for every element of a N-D tensor array. More...
 
template<... >
void tiffToRGBImage (...)
 Import a RGB (3-band or color-mapped) TiffImage into a RGB image. More...
 
template<... >
void tiffToScalarImage (...)
 Convert single-band TiffImage to scalar image. More...
 
std::string tolower (std::string s)
 
template<... >
void totalVariationFilter (...)
 Performs standard Total Variation Regularization. More...
 
template<... >
void transformImage (...)
 Apply unary point transformation to each pixel. More...
 
template<... >
void transformImageIf (...)
 Apply unary point transformation to each pixel within the ROI (i.e., where the mask is non-zero). More...
 
template<... >
void transformMultiArray (...)
 Transform a multi-dimensional array with a unary function or functor. More...
 
template<class SplineImage , class DestIterator , class DestAccessor >
void transformToPolarCoordinates (SplineImage const &src, DestIterator d_ul, DestIterator d_lr, DestAccessor d_acc)
 Transforms a given image to its (image-centered) polar coordinates representation. More...
 
linalg::TemporaryMatrix< double > translationMatrix2D (TinyVector< double, 2 > const &shift)
 Create homogeneous matrix representing a 2D translation. More...
 
template<class V , int SIZE, class T >
TinyVector< V, SIZE > transpose (TinyVector< V, SIZE > const &t, TinyVector< T, SIZE > const &permutation)
 transposed copy More...
 
template<... >
void transposeImage (...)
 Transpose an image over the major or minor diagonal. More...
 
void uncompress (char const *source, std::size_t srcSize, char *dest, std::size_t destSize, CompressionMethod method)
 
template<... >
unsigned int unionFindWatershedsBlockwise (...)
 Blockwise union-find watersheds transform for MultiArrays and ChunkedArrays. More...
 
template<... >
void vectorToTensor (...)
 Calculate the tensor (outer) product of a 2D vector with itself. More...
 
template<... >
void vectorToTensorMultiArray (...)
 Calculate the tensor (outer) product of a N-D vector with itself. More...
 
template<... >
unsigned int watersheds3D (...)
 Region Segmentation by means of the watershed algorithm. More...
 
template<... >
Label watershedsMultiArray (...)
 Watershed segmentation of an arbitrary-dimensional array. More...
 
template<... >
unsigned int watershedsRegionGrowing (...)
 Region segmentation by means of a flooding-based watershed algorithm. More...
 
template<... >
unsigned int watershedsUnionFind (...)
 Region segmentation by means of the union-find watershed algorithm. More...
 
template<... >
void writeHDF5 (...)
 Store array data in an HDF5 file. More...
 
template<size_t N, class T , class C >
void writeHDF5Attr (hid_t loc, const char *name, MultiArrayView< N, T, C > const &array)
 
template<size_t N, class C >
void writeHDF5Attr (hid_t loc, const char *name, MultiArrayView< N, std::string, C > const &array)
 
template<class T >
void writeHDF5Attr (hid_t loc, const char *name, ArrayVectorView< T > &array)
 
template<class Arr >
void writeHDF5Attr (std::string filePath, std::string pathInFile, Arr &ar)
 
template<class V >
TinyVector< float, 3 > yPrimeCbCr2Polar (V const &ycbcr)
 Create polar representation form Y'CbCr. More...
 
template<class V >
TinyVector< float, 3 > yPrimeIQ2Polar (V const &yiq)
 Create polar representation form Y'IQ. More...
 
template<class V >
TinyVector< float, 3 > yPrimePbPr2Polar (V const &ypbpr)
 Create polar representation form Y'PbPr. More...
 
template<class V >
TinyVector< float, 3 > yPrimeUV2Polar (V const &yuv)
 Create polar representation form Y'UV. More...
 
template<... >
void applyWindowFunction (...)
 Apply a window function to each pixels of a given image. More...
 
template<... >
void medianFilter (...)
 This function calculates the median of a window of given size for the complete image. More...
 
template<... >
void upwindImage (...)
 This function calculates discrete upwinding scheme proposed by J. Weickert (2002) in Coherence-Enhancing Show Filters. More...
 
template<... >
void frostFilter (...)
 This function tries to reduce the speckle noise of an image by applying the basic Frost filter. More...
 

Detailed Description

all VIGRA functionality is located in namespace vigra

This header provides definitions of graph-related algorithms

This header provides definitions of graph-related types and optionally provides a gateway to popular graph libraries (for now, BGL is supported).

Typedef Documentation

Standard ridge regression split

Enumeration Type Documentation

tags used with the RandomForestOptions class

See Also
RF_Traits::Option_t

Function Documentation

void vigra::applyWindowFunction (   ...)

Apply a window function to each pixels of a given image.

This function calculates the results for a window function (given as a functor) when applied to the complete image. Also allows a correct border handling! See medianFilter() for an example of a quite basic window function and its application.If you pass a functor to this function, which implements the two functions:

  1. Diff2D windowShape() const
    to return the filter window size, which has to be odd in each dimension and
  2. void operator()(SrcIterator s, SrcAccessor s_acc, DestIterator d, DestAccessor d_acc)
    to compute the results of the current window,

this function calculates the results for the complete image.

All border treatment modes (except BORDER_TREATMENT_CLIP) are supported.

The input pixel type T1 must be a linear space over the window functions' value_type T, i.e. addition of source values, multiplication with functions' values, and NumericTraits must be defined. The filters' value_type must be an algebraic field, i.e. the arithmetic operations (+, -, *, /) and NumericTraits must be defined.

Declarations:

pass 2D array views:

namespace vigra {
template <class T1, class S1,
class T2, class S2,
class ProcessingFunctor>
void
applyWindowFunction(MultiArrayView<2, T1, S1> const & src,
MultiArrayView<2, T2, S2> dest,
ProcessingFunctor func,
BorderTreatmentMode border = BORDER_TREATMENT_REPEAT);
}

show deprecated declarations

Usage:

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

unsigned int w=1000, h=1000;
MultiArray<2, float> src(w,h), dest(w,h);
...
template<class VALUETYPE>
class AvgFunctor
{
public:
MedianFunctor(Diff2D window_shape)
: m_window_shape(window_shape)
{
}
template <class SrcIterator, class SrcAccessor, class DestIterator, class DestAccessor>
void operator()(SrcIterator s, SrcAccessor s_acc, DestIterator d, DestAccessor d_acc)
{
SrcIterator s_ul = s - m_window_shape/2,
s_lr = s_ul + m_window_shape;
VALUETYPE result = NumericTraits<float>::zero();
SrcIterator ys = s_ul;
SrcIterator xs = ys;
for( ; ys.y != s_lr.y; ys.y++)
{
for(xs = ys; xs.x != s_lr.x; xs.x++, iter++)
{
res += s_acc(xs);
}
}
d_acc.set(res/(m_window_shape.x*m_window_shape.y),d);
}
Diff2D windowShape() const
{
return m_window_shape;
}
private:
Diff2D m_window_shape;
};
// create an AverageFilter function for a 5x5 filter
AvgFunctor func(Diff2D(5,5));
// apply the filter function to the input image
applyWindowFunction(src, dest, func);

Preconditions:

The image must be larger than the window size.

vigra::doxygen_overloaded_function ( template<...> void  separableConvolveBlockwise)

Separated convolution on ChunkedArrays.

Declarations:

namespace vigra {
// apply each kernel from the sequence 'kernels' in turn
template <unsigned int N, class T1, class T2, class KernelIterator>
void separableConvolveBlockwise(const ChunkedArra<N, T1>& source, ChunkedArray<N, T2>& destination, KernelIterator kernels);
// apply the same kernel to all dimensions
template <unsigned int N, class T1, class T2, class T3>
void separableConvolveBlockwise(const ChunkedArra<N, T1>& source, ChunkedArray<N, T2>& destination, Kernel1D<T3> const & kernel);
}

This function computes a separated convolution for a given ChunkedArray. For infinite precision T1, this is equivalent to separableConvolveMultiArray. In practice, floating point inaccuracies will make the result differ slightly.

void vigra::compress ( char const *  source,
std::size_t  size,
ArrayVector< char > &  dest,
CompressionMethod  method 
)

Compress the source buffer.

The destination array will be resized as required.

void vigra::uncompress ( char const *  source,
std::size_t  srcSize,
char *  dest,
std::size_t  destSize,
CompressionMethod  method 
)

Uncompress the source buffer when the uncompressed size is known.

The destination buffer must be allocated to the correct size.

void vigra::projectBack ( const AdjacencyListGraph &  rag,
const BASE_GRAPH &  bg,
const Int64  ignoreLabel,
const BASE_GRAPH_LABELS  bgLabels,
const RAG_FEATURES &  ragFeatures,
BASE_GRAPH_FEATURES &  bgFeatures 
)

project node features of a region adjacency graph back to the base graph.

This function can be used to show a segmentation or node features of RAG on pixel / voxel level

void vigra::medianFilter (   ...)

This function calculates the median of a window of given size for the complete image.

This function calculates the median of a window of given size for the complete image. It also allows a correct border handling, since it uses the applyWindowFunction environment for computation!All border treatment modes (except BORDER_TREATMENT_CLIP) are supported.

The input pixel type T1 must be a linear space over the window functions' value_type T. Especially, the values must be sortable by std::sort, to derive the mean values aka the median.

Declarations:

pass 2D array views:

namespace vigra {
template <class T1, class S1,
class T2, class S2>
void
medianFilter(MultiArrayView<2, T1, S1> const & src,
MultiArrayView<2, T2, S2> dest,
Diff2D window_shape,
BorderTreatmentMode border = BORDER_TREATMENT_REPEAT);
}

show deprecated declarations

Usage:

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

unsigned int w=1000, h=1000;
MultiArray<2, float> src(w,h), dest(w,h);
...
// apply a median filter with a window size of 5x5
medianFilter(src, dest, Diff2D(5,5));

Preconditions:

The image must be larger than the window size of the filter.

detail::RF_DEFAULT & rf_default ( )

factory function to return a RF_DEFAULT tag

See Also
RandomForest<>::learn()
void vigra::rf_export_HDF5 ( const RandomForest< T, Tag > &  rf,
HDF5File &  h5context,
const std::string &  pathname = "" 
)

Save a random forest to an HDF5File object into a specified HDF5 group.

The random forest is saved as a set of HDF5 datasets, groups, and attributes below a certain HDF5 group (default: current group of the HDF5File object). No additional data should be stored in that group.

Parameters
rfRandom forest object to be exported
h5contextHDF5File object to use
pathnameIf empty or not supplied, save the random forest to the current group of the HDF5File object. Otherwise, save to a new-created group specified by the path name, which may be either relative or absolute.
void vigra::rf_export_HDF5 ( const RandomForest< T, Tag > &  rf,
const std::string &  filename,
const std::string &  pathname = "" 
)

Save a random forest to a named HDF5 file into a specified HDF5 group.

The random forest is saved as a set of HDF5 datasets, groups, and attributes below a certain HDF5 group (default: root). No additional data should be stored in that group.

Parameters
rfRandom forest object to be exported
filenameName of an HDF5 file to open
pathnameIf empty or not supplied, save the random forest to the root group of the HDF5 file. Otherwise, save to a new-created group specified by the path name (relative to the root group).
void vigra::rf_export_HDF5 ( const RandomForest< T, Tag > &  rf,
hid_t  outf_id,
const std::string &  pathname = "" 
)

Save a random forest to an HDF5 file specified by its id.

The random forest is saved as a set of HDF5 datasets, groups, and attributes below a certain HDF5 group (default: root). No additional data should be stored in that group.

Warning
In case the underlying HDF5 library used by Vigra is not exactly the same library used to open the file with the given id, this method will lead to crashes.
Parameters
rfRandom forest object to be exported
outf_idHDF5 file id
pathnameIf empty or not supplied, save the random forest to the root group of the HDF5 file. Otherwise, save to a new-created group specified by the path name (relative to the root group).
bool vigra::rf_import_HDF5 ( RandomForest< T, Tag > &  rf,
HDF5File &  h5context,
const std::string &  pathname = "" 
)

Read a random forest from an HDF5File object's specified group.

The random forest is read from a certain HDF5 group (default: current group of the HDF5File object) as a set of HDF5 datasets, groups, and attributes. No additional data should be present in that group.

Parameters
rfRandom forest object to be imported
h5contextHDF5File object to use
pathnameIf empty or not supplied, read from the random forest from the current group of the HDF5File object. Otherwise, use the group specified by the path name, which may be either relative or absolute.
bool vigra::rf_import_HDF5 ( RandomForest< T, Tag > &  rf,
const std::string &  filename,
const std::string &  pathname = "" 
)

Read a random forest from a named HDF5 file's specified group.

The random forest is read from a certain HDF5 group (default: root group of the HDF5 file) as a set of HDF5 datasets, groups, and attributes. No additional data should be present in that group.

Parameters
rfRandom forest object to be imported
filenameName of an HDF5 file to open
pathnameIf empty or not supplied, read from the random forest from the current group of the HDF5 file. Otherwise, use the group specified by the path name, which may be either relative or absolute.
bool vigra::rf_import_HDF5 ( RandomForest< T, Tag > &  rf,
hid_t  inf_id,
const std::string &  pathname = "" 
)

Read a random forest from an HDF5 file specified by its id.

The random forest is read from a certain HDF5 group (default: root group of the HDF5 file) as a set of HDF5 datasets, groups, and attributes. No additional data should be present in that group.

Warning
In case the underlying HDF5 library used by Vigra is not exactly the same library used to open the file with the given id, this method will lead to crashes.
Parameters
rfRandom forest object to be imported
inf_idHDF5 file id
pathnameIf empty or not supplied, read from the random forest from the current group of the HDF5 file. Otherwise, use the group specified by the path name, which may be either relative or absolute.
void upwindImage (   ...)

This function calculates discrete upwinding scheme proposed by J. Weickert (2002) in Coherence-Enhancing Show Filters.

This function calculates of the coherence enhancing shock filter proposed by J. Weickert (2002) in Coherence-Enhancing Show Filters.

This function calculates of the coherence enhancing shock filter proposed by J. Weickert (2002): Coherence-Enhancing Show Filters. It uses the structure tensor information of an image and an iterative discrete upwinding scheme instead of pure dilation and erosion to perform the necessary morphological operations on the image.An image is upwinded positively (dilation), if the given second src is positive. Otherwise it is upwinds negatively (eroded). The effect can be steered by an upwinding factor.

Declarations:

pass 2D array views:

namespace vigra {
template <class T1, class S1,
class T2, class S2
class T3, class S3>
void
upwindImage(MultiArrayView<2, T1, S1> const & src,
MultiArrayView<2, T2, S2> const & src2,
MultiArrayView<2, T3, S3> dest,
float upwind_factor_h);
}

show deprecated declarations

Declarations:

pass 2D array views:

namespace vigra {
template <class T1, class S1,
class T2, class S2>
void
shockFilter(MultiArrayView<2, T1, S1> const & src,
MultiArrayView<2, T2, S2> dest,
float sigma, float rho, float upwind_factor_h,
unsigned int iterations);
}

show deprecated declarations

Usage:

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

unsigned int w=1000, h=1000;
MultiArray<2, float> src(w,h), dest(w,h);
...
// apply a shock-filter:
shockFilter(src, dest, 1.0, 5.0, 0.3, 5);

Preconditions:

The image must be larger than the window size of the filter.

void vigra::frostFilter (   ...)

This function tries to reduce the speckle noise of an image by applying the basic Frost filter.

This function tries to reduce the speckle noise of an image by means of applying the basic Frost filter using a window of given size and a damping factor k. The implementation is according to the article by Lopez & Touzi & Nezry (1990): Adaptive speckle filters and scene heterogenity.The user has to provide a window size and a damping factor k. The implementation is according to the article by Lopez & Touzi & Nezry (1990): Adaptive speckle filters and scene heterogenity.

All restrictions of the called functions applyWindowFunction apply.

Preconditions:

0.0 < k <= 1.0

Declarations:

pass 2D array views:

namespace vigra {
template <class T1, class S1,
class T2, class S2>
void
frostFilter(MultiArrayView<2, T1, S1> const & src,
MultiArrayView<2, T2, S2> dest,
Diff2D window_shape, float k,
BorderTreatmentMode border = BORDER_TREATMENT_REPEAT);
}

show deprecated declarations

Usage:

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

unsigned int w=1000, h=1000;
MultiArray<2, float> src(w,h), dest(w,h);
...
// apply a (basic) frost filter with a window size of 5x5 and a damping factor of 0.5
frostFilter(src, dest, Diff2D(5,5), 0.5);
std::string vigra::asString ( t)

Convert a value to a string. Available for integral and floating point types and void *.

std::string vigra::tolower ( std::string  s)

Convert string to lower case.

std::string vigra::normalizeString ( std::string const &  s)

Convert string to lower case and remove any white space characters.

void vigra::for_each_in_tuple ( TPL &&  t,
FUNCTOR &&  f 
)

The for_each_in_tuple function calls the functor f on all elements of the tuple t. For each element type in the tuple, the functor must have an appropriate overload of operator().

Example:

#include <iostream>
#include <tuple>
struct print
{
template <typename T>
void operator()(T const & t) const
{
std::cout << t << std::endl;
}
};
struct add_one
{
template <typename T>
void operator()(T & t) const
{
t += 1;
}
};
int main()
{
std::tuple<int, double, size_t> tpl(-5, 0.4, 10);
vigra::for_each_in_tuple(tpl, add_one());
vigra::for_each_in_tuple(tpl, print());
}

© 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)