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

details Mathematical Functions VIGRA

Namespaces

 vigra
 

Classes

class  BSpline< ORDER, T >
 
class  BSplineBase< ORDER, T >
 
class  CatmullRomSpline< T >
 
class  CoscotFunction< T >
 
class  Gaussian< T >
 

Functions

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 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...
 
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...
 
UInt32 ceilPower2 (UInt32 x)
 Round up to the nearest power of 2. More...
 
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 T1 , class T2 >
bool closeAtTolerance (T1 l, T2 r, typename PromoteTraits< T1, T2 >::Promote epsilon)
 Tolerance based floating-point equality. More...
 
UInt32 concatenateChecksum (UInt32 checksum, const char *data, unsigned int size)
 
template<class REAL >
REAL cos_pi (REAL x)
 cos(pi*x). 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...
 
bool even (int t)
 Check if an integer is even. More...
 
UInt32 floorPower2 (UInt32 x)
 Round down to the nearest power of 2. More...
 
double gamma (double x)
 The gamma function. More...
 
template<typename IntType >
IntType gcd (IntType n, IntType m)
 
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<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 inspectSequence (...)
 Call an analyzing functor at every element of a sequence. More...
 
template<class InIterator , class OutIterator >
void inversePermutation (InIterator first, InIterator last, OutIterator out)
 Compute the inverse of a given permutation. 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<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 Iterator , class Value >
void linearSequence (Iterator first, Iterator last, Value start, Value step)
 Fill an array with a sequence of numbers. 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 ArrayLike , class Compare >
detail::IndexCompare
< ArrayLike, Compare > 
makeIndexComparator (ArrayLike a, Compare c)
 Create a compare functor for indirect sort. 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<class T >
NormTraits< T >::NormType norm (T const &t)
 The norm of a numerical object. More...
 
bool odd (int t)
 Check if an integer is odd. More...
 
result_type operator() (argument_type x) const
 
template<unsigned n, class V >
power (const V &x)
 Exponentiation to a positive integer power by squaring. More...
 
REAL round (REAL v)
 The rounding function. More...
 
long long roundi (double t)
 Round and cast to integer. More...
 
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<class REAL >
REAL sin_pi (REAL x)
 sin(pi*x). More...
 
template<class T >
NumericTraits< T >::Promote sq (T t)
 The square function. More...
 
Int32 sqrti (Int32 v)
 Signed integer square root. More...
 
NormTraits< T >::SquaredNormType squaredNorm (T const &t)
 The squared norm of a numerical object. More...
 
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...
 

Detailed Description

Useful mathematical functions and functors.

Function Documentation

Iterator vigra::argMin ( Iterator  first,
Iterator  last 
)

Find the minimum element in a sequence.

The function returns the iterator referring to the minimum element. This is identical to the function std::min_element().

Required Interface:

Iterator is a standard forward iterator.
bool f = *first < NumericTraits<typename std::iterator_traits<Iterator>::value_type>::max();

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

Iterator vigra::argMax ( Iterator  first,
Iterator  last 
)

Find the maximum element in a sequence.

The function returns the iterator referring to the maximum element. This is identical to the function std::max_element().

Required Interface:

Iterator is a standard forward iterator.
bool f = NumericTraits<typename std::iterator_traits<Iterator>::value_type>::min() < *first;

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

Iterator vigra::argMinIf ( Iterator  first,
Iterator  last,
UnaryFunctor  condition 
)

Find the minimum element in a sequence conforming to a condition.

The function returns the iterator referring to the minimum element, where only elements conforming to the condition (i.e. where condition(*iterator) evaluates to true) are considered. If no element conforms to the condition, or the sequence is empty, the end iterator last is returned.

Required Interface:

Iterator is a standard forward iterator.
bool c = condition(*first);
bool f = *first < NumericTraits<typename std::iterator_traits<Iterator>::value_type>::max();

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

Iterator vigra::argMaxIf ( Iterator  first,
Iterator  last,
UnaryFunctor  condition 
)

Find the maximum element in a sequence conforming to a condition.

The function returns the iterator referring to the maximum element, where only elements conforming to the condition (i.e. where condition(*iterator) evaluates to true) are considered. If no element conforms to the condition, or the sequence is empty, the end iterator last is returned.

Required Interface:

Iterator is a standard forward iterator.
bool c = condition(*first);
bool f = NumericTraits<typename std::iterator_traits<Iterator>::value_type>::min() < *first;

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

void vigra::linearSequence ( Iterator  first,
Iterator  last,
Value  start,
Value  step 
)

Fill an array with a sequence of numbers.

The sequence starts at start and is incremented with step. Default start and stepsize are 0 and 1 respectively. This is a generalization of function std::iota() in C++11.

Declaration:

namespace vigra {
template <class Iterator, class Value>
void linearSequence(Iterator first, Iterator last,
Value const & start = 0, Value const & step = 1);
}

Required Interface:

Iterator is a standard forward iterator.
first = start;
start += step;

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

void vigra::inspectSequence (   ...)

Call an analyzing functor at every element of a sequence.

This function can be used to collect statistics of the sequence [first, last) defined by these two input interators. The results must be stored in the functor, which serves as a return value.

Declarations:

namespace vigra {
template <class InputIterator, class Functor>
void
inspectSequence(InputIterator first, InputIterator last, Functor & f);
}

Usage:

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

std::vector array(100);
// init functor
vigra::inspectSequence(array.begin(), array.end(), minmax);
cout << "Min: " << minmax.min << " Max: " << minmax.max;
detail::IndexCompare<ArrayLike, Compare> vigra::makeIndexComparator ( ArrayLike  a,
Compare  c 
)

Create a compare functor for indirect sort.

Indirect sorting refers to the situation where you have an array holding data and another array holding indices referencing the first array, and you want to sort the index array according to some property of the data array without changing the data array itself. The factory function makeIndexComparator() creates a sorting predicate for this task, given a sorting predicate for the data array.

See Also
vigra::indexSort(), vigra::applyPermutation()

Usage:

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

const std:vector<double> data(...); // data is immutable
std::vector<int> index(data.size());
std::iota(index.begin(), index.end());
// sort the indices such that data[index[k]] is an ascending sequence in k
std::sort(index.begin(), index.end(), makeIndexComparator(data));

Declarations:

namespace vigra {
// compare using std::less
template <class ArrayLike>
auto makeIndexComparator(ArrayLike a);
// compare using functor Compare
template <class ArrayLike, class Compare>
auto makeIndexComparator(ArrayLike a, Compare c);
}
void vigra::indexSort ( Iterator  first,
Iterator  last,
IndexIterator  index_first,
Compare  c 
)

Return the index permutation that would sort the input array.

To actually sort an array according to the ordering thus determined, use applyPermutation().

Usage:

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

const std:vector<double> data(...); // data is immutable
std::vector<int> index(data.size());
// arrange indices such that data[index[k]] is an ascending sequence in k
indexSort(data.begin(), data.end(), index.begin());

Declarations:

namespace vigra {
// compare using std::less
template <class Iterator, class IndexIterator>
void indexSort(Iterator first, Iterator last, IndexIterator index_first);
// compare using functor Compare
template <class Iterator, class IndexIterator, class Compare>
void indexSort(Iterator first, Iterator last, IndexIterator index_first, Compare compare);
}

Required Interface:

Iterator and IndexIterators are random access iterators.
bool res = compare(first[*index_first], first[*index_first]);

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

void vigra::applyPermutation ( IndexIterator  index_first,
IndexIterator  index_last,
InIterator  in,
OutIterator  out 
)

Sort an array according to the given index permutation.

The iterators in and out may not refer to the same array, as this would overwrite the input prematurely.

Declaration:

namespace vigra {
template <class IndexIterator, class InIterator, class OutIterator>
void applyPermutation(IndexIterator index_first, IndexIterator index_last,
InIterator in, OutIterator out);
}

Required Interface:

OutIterator and IndexIterators are forward iterators.
InIterator is a random access iterator.
out = in[*index_first];

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

void vigra::inversePermutation ( InIterator  first,
InIterator  last,
OutIterator  out 
)

Compute the inverse of a given permutation.

This is just another name for indexSort(), referring to another semantics.

Declaration:

namespace vigra {
template <class InIterator, class OutIterator>
void inversePermutation(InIterator first, InIterator last,
OutIterator out);
}

Required Interface:

InIterator and OutIterator are random access iterators.
out = in[*index_first];

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

UInt32 vigra::checksum ( const char *  data,
unsigned int  size 
)

Compute the CRC-32 checksum of a byte array.

Implementation note: This function is slower on big-endian machines because the "4 bytes at a time" optimization is only implemented for little-endian.

UInt32 vigra::concatenateChecksum ( UInt32  checksum,
const char *  data,
unsigned int  size 
)

Concatenate a byte array to an existing CRC-32 checksum.

double vigra::besselJ ( int  n,
double  x 
)

Bessel function of the first kind.

Computes the value of BesselJ of integer order n and argument x. Negative x are unsupported and will result in a std::domain_error.

This function wraps a number of existing implementations and falls back to a rather slow algorithm if none of them is available. In particular, it uses boost::math when HasBoostMath is #defined, or native implementations on gcc and MSVC otherwise.

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

double vigra::besselY ( int  n,
double  x 
)

Bessel function of the second kind.

Computes the value of BesselY of integer order n and argument x. Negative x are unsupported and will result in a std::domain_error.

This function wraps a number of existing implementations and falls back to a rather slow algorithm if none of them is available. In particular, it uses boost::math when HasBoostMath is #defined, or native implementations on gcc and MSVC otherwise.

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

REAL vigra::round ( REAL  v)

The rounding function.

Defined for all floating point types. Rounds towards the nearest integer such that abs(round(t)) == round(abs(t)) for all t.

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

long long vigra::roundi ( double  t)

Round and cast to integer.

Rounds to the nearest integer like round(), but casts the result to long long (this will be faster and is usually needed anyway).

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

UInt32 vigra::ceilPower2 ( UInt32  x)

Round up to the nearest power of 2.

Efficient algorithm for finding the smallest power of 2 which is not smaller than x (function clp2() from Henry Warren: "Hacker's Delight", Addison-Wesley, 2003, see http://www.hackersdelight.org/). If x > 2^31, the function will return 0 because integer arithmetic is defined modulo 2^32.

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

UInt32 vigra::floorPower2 ( UInt32  x)

Round down to the nearest power of 2.

Efficient algorithm for finding the largest power of 2 which is not greater than x (function flp2() from Henry Warren: "Hacker's Delight", Addison-Wesley, 2003, see http://www.hackersdelight.org/).

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

Int32 vigra::log2i ( UInt32  x)

Compute the base-2 logarithm of an integer.

Returns the position of the left-most 1-bit in the given number x, or -1 if x == 0. That is,

assert(k >= 0 && k < 32 && log2i(1 << k) == k);

The function uses Robert Harley's algorithm to determine the number of leading zeros in x (algorithm nlz10() at http://www.hackersdelight.org/). But note that the functions floorPower2() or ceilPower2() are more efficient and should be preferred when possible.

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

NumericTraits<T>::Promote vigra::sq ( t)

The square function.

sq(x) = x*x is needed so often that it makes sense to define it as a function.

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

V vigra::power ( const V &  x)

Exponentiation to a positive integer power by squaring.

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

UInt32 sqrti ( Int32  v)

Signed integer square root.

Unsigned integer square root.

Useful for fast fixed-point computations.

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

T vigra::sign ( t)

The sign function.

Returns 1, 0, or -1 depending on the sign of t, but with the same type as t.

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

int vigra::signi ( t)

The integer sign function.

Returns 1, 0, or -1 depending on the sign of t, converted to int.

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

T1 vigra::sign ( T1  t1,
T2  t2 
)

The binary sign function.

Transfers the sign of t2 to t1.

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

bool vigra::even ( int  t)

Check if an integer is even.

Defined for all integral types.

bool vigra::odd ( int  t)

Check if an integer is odd.

Defined for all integral types.

NormTraits<T>::SquaredNormType vigra::squaredNorm ( T const &  t)

The squared norm of a numerical object.

  • For scalar types: equals vigra::sq(t).
  • For vectorial types (including TinyVector): equals vigra::dot(t, t).
  • For complex number types: equals vigra::sq(t.real()) + vigra::sq(t.imag()).
  • For array and matrix types: results in the squared Frobenius norm (sum of squares of the matrix elements).
NormTraits<T>::NormType vigra::norm ( T const &  t)

The norm of a numerical object.

For scalar types: implemented as abs(t)
otherwise: implemented as sqrt(squaredNorm(t)).

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

void vigra::symmetric2x2Eigenvalues ( a00,
a01,
a11,
T *  r0,
T *  r1 
)

Compute the eigenvalues of a 2x2 real symmetric matrix.

This uses the analytical eigenvalue formula

\[ \lambda_{1,2} = \frac{1}{2}\left(a_{00} + a_{11} \pm \sqrt{(a_{00} - a_{11})^2 + 4 a_{01}^2}\right) \]

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

void vigra::symmetric3x3Eigenvalues ( a00,
a01,
a02,
a11,
a12,
a22,
T *  r0,
T *  r1,
T *  r2 
)

Compute the eigenvalues of a 3x3 real symmetric matrix.

This uses a numerically stable version of the analytical eigenvalue formula according to

David Eberly: "Eigensystems for 3 x 3 Symmetric Matrices (Revisited)", Geometric Tools Documentation, 2006

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

double vigra::ellipticIntegralF ( double  x,
double  k 
)

The incomplete elliptic integral of the first kind.

This function computes

\[ \mbox{F}(x, k) = \int_0^x \frac{1}{\sqrt{1 - k^2 \sin(t)^2}} dt \]

according to the algorithm given in Press et al. "Numerical Recipes".

Note: In some libraries (e.g. Mathematica), the second parameter of the elliptic integral functions must be k^2 rather than k. Check the documentation when results disagree!

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

double vigra::ellipticIntegralE ( double  x,
double  k 
)

The incomplete elliptic integral of the second kind.

This function computes

\[ \mbox{E}(x, k) = \int_0^x \sqrt{1 - k^2 \sin(t)^2} dt \]

according to the algorithm given in Press et al. "Numerical Recipes". The complete elliptic integral of the second kind is simply ellipticIntegralE(M_PI/2, k).

Note: In some libraries (e.g. Mathematica), the second parameter of the elliptic integral functions must be k^2 rather than k. Check the documentation when results disagree!

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

double vigra::chi2 ( unsigned int  degreesOfFreedom,
double  arg,
double  accuracy = 1e-7 
)

Chi square distribution.

Computes the density of a chi square distribution with degreesOfFreedom and tolerance accuracy at the given argument arg by calling noncentralChi2(degreesOfFreedom, 0.0, arg, accuracy).

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

double vigra::chi2CDF ( unsigned int  degreesOfFreedom,
double  arg,
double  accuracy = 1e-7 
)

Cumulative chi square distribution.

Computes the cumulative density of a chi square distribution with degreesOfFreedom and tolerance accuracy at the given argument arg, i.e. the probability that a random number drawn from the distribution is below arg by calling noncentralChi2CDF(degreesOfFreedom, 0.0, arg, accuracy).

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

double vigra::noncentralChi2 ( unsigned int  degreesOfFreedom,
double  noncentrality,
double  arg,
double  accuracy = 1e-7 
)

Non-central chi square distribution.

Computes the density of a non-central chi square distribution with degreesOfFreedom, noncentrality parameter noncentrality and tolerance accuracy at the given argument arg. It uses Algorithm AS 231 from Appl. Statist. (1987) Vol.36, No.3 (code ported from http://lib.stat.cmu.edu/apstat/231). The algorithm has linear complexity in the number of degrees of freedom.

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

double vigra::noncentralChi2CDF ( unsigned int  degreesOfFreedom,
double  noncentrality,
double  arg,
double  accuracy = 1e-7 
)

Cumulative non-central chi square distribution.

Computes the cumulative density of a chi square distribution with degreesOfFreedom, noncentrality parameter noncentrality and tolerance accuracy at the given argument arg, i.e. the probability that a random number drawn from the distribution is below arg It uses Algorithm AS 231 from Appl. Statist. (1987) Vol.36, No.3 (code ported from http://lib.stat.cmu.edu/apstat/231). The algorithm has linear complexity in the number of degrees of freedom (see noncentralChi2CDFApprox() for a constant-time algorithm).

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

double vigra::noncentralChi2CDFApprox ( unsigned int  degreesOfFreedom,
double  noncentrality,
double  arg 
)

Cumulative non-central chi square distribution (approximate).

Computes approximate values of the cumulative density of a chi square distribution with degreesOfFreedom, and noncentrality parameter noncentrality at the given argument arg, i.e. the probability that a random number drawn from the distribution is below arg It uses the approximate transform into a normal distribution due to Wilson and Hilferty (see Abramovitz, Stegun: "Handbook of Mathematical Functions", formula 26.3.32). The algorithm's running time is independent of the inputs, i.e. is should be used when noncentralChi2CDF() is too slow, and approximate values are sufficient. The accuracy is only about 0.1 for few degrees of freedom, but reaches about 0.001 above dof = 5.

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

REAL vigra::legendre ( unsigned int  l,
int  m,
REAL  x 
)

Associated Legendre polynomial.

Computes the value of the associated Legendre polynomial of order l, m for argument x. x must be in the range [-1.0, 1.0], otherwise an exception is thrown. The standard Legendre polynomials are the special case m == 0.

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

REAL vigra::legendre ( unsigned int  l,
REAL  x 
)

Legendre polynomial.

Computes the value of the Legendre polynomial of order l for argument x. x must be in the range [-1.0, 1.0], otherwise an exception is thrown.

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

REAL vigra::sin_pi ( REAL  x)

sin(pi*x).

Essentially calls std::sin(M_PI*x) but uses a more accurate implementation to make sure that sin_pi(1.0) == 0.0 (which does not hold for std::sin(M_PI) due to round-off error), and sin_pi(0.5) == 1.0.

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

REAL vigra::cos_pi ( REAL  x)

cos(pi*x).

Essentially calls std::cos(M_PI*x) but uses a more accurate implementation to make sure that cos_pi(1.0) == -1.0 and cos_pi(0.5) == 0.0.

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

double vigra::gamma ( double  x)

The gamma function.

This function implements the algorithm from
Zhang and Jin: "Computation of Special Functions", John Wiley and Sons, 1996.

The argument must be <= 171.0 and cannot be zero or a negative integer. An exception is thrown when these conditions are violated.

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

Examples:
total_variation.cxx.
double vigra::loggamma ( double  x)

The natural logarithm of the gamma function.

This function is based on a free implementation by Sun Microsystems, Inc., see sourceware.org archive. It can be removed once all compilers support the new C99 math functions.

The argument must be positive and < 1e30. An exception is thrown when these conditions are violated.

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

bool vigra::closeAtTolerance ( T1  l,
T2  r,
typename PromoteTraits< T1, T2 >::Promote  epsilon 
)

Tolerance based floating-point equality.

Check whether two floating point numbers are equal within the given tolerance. This is useful because floating point numbers that should be equal in theory are rarely exactly equal in practice. If the tolerance epsilon is not given, twice the machine epsilon is used.

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

bool vigra::lessEqualAtTolerance ( T1  l,
T2  r,
typename PromoteTraits< T1, T2 >::Promote  epsilon 
)

Tolerance based floating-point less-or-equal.

Check whether two floating point numbers are less or equal within the given tolerance. That is, l can actually be greater than r within the given epsilon. This is useful because floating point numbers that should be equal in theory are rarely exactly equal in practice. If the tolerance epsilon is not given, twice the machine epsilon is used.

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

bool vigra::greaterEqualAtTolerance ( T1  l,
T2  r,
typename PromoteTraits< T1, T2 >::Promote  epsilon 
)

Tolerance based floating-point greater-or-equal.

Check whether two floating point numbers are greater or equal within the given tolerance. That is, l can actually be less than r within the given epsilon. This is useful because floating point numbers that should be equal in theory are rarely exactly equal in practice. If the tolerance epsilon is not given, twice the machine epsilon is used.

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

IntType vigra::gcd ( IntType  n,
IntType  m 
)

Calculate the greatest common divisor.

This function works for arbitrary integer types, including user-defined (e.g. infinite precision) ones.

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

IntType vigra::lcm ( IntType  n,
IntType  m 
)

Calculate the lowest common multiple.

This function works for arbitrary integer types, including user-defined (e.g. infinite precision) ones.

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

CatmullRomSpline< T >::result_type operator() ( argument_type  x) const

function (functor) call

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