'utils' - a collection of useful math and programming tools


Classes

class  Go::Array< T, Dim >
 Compile-time sized array. More...
class  Go::BaryCoordSystem< Dim >
 Encapsulates a barycentric coordinatesystem. More...
class  Go::BaryCoordSystemTriangle3D
 A barycentric coordinate system for a triangle (2-manifold) embedded in 3D. More...
class  Go::BoundingBox
 Axis-aligned bounding box. More...
class  Go::Fun2Fun< Functor >
 Brief description. More...
class  Go::CompositeBox
 Composite of two bounding boxes. More...
class  Go::CoordinateSystem< Dim >
 Defines a Cartesian coordinate system. More...
class  Go::CPUclock
 A class for measuring CPU time in programs. More...
class  Go::DirectionCone
 Class representing a direction cone in Euclidian N-space The direction cone is characterised by a central direction and an angle 'alpha'. More...
struct  Go::Factorial< N >
 Calculate the factorial of a given integer N. More...
struct  Go::Factorial< 1 >
struct  Go::InverseFactorial< T, N >
 Compute the inverse of the factorial of a given integer N, where T is typically 'float' or 'double'. More...
class  Go::FunctionMinimizer< Functor >
 This is the FunctionMinimizer class that can be used ex. More...
class  Go::MatrixXD< T, Dim >
 Square n-dimensional (compile time constant dim) matrix. More...
class  Go::Point
 Run-time sized point class. More...
class  Go::Rational
 Class representing rational numbers. More...
class  Go::RotatedBox
 A rotated version of CompositeBox. More...
class  Go::ScratchVect< T, N >
 A template Vector class that behaves much like the std::vector template, but stores its elements on the stack rather than on the heap as long as the total size stays less than 'N' (template argument). More...

Typedefs

typedef Array< double, 2 > Go::Vector2D
 Typedef for ease of use in frequently used case.
typedef Array< double, 3 > Go::Vector3D
 Typedef for ease of use in frequently used case.
typedef Array< double, 4 > Go::Vector4D
 Typedef for ease of use in frequently used case.
typedef BaryCoordSystem< 2 > Go::BaryCoordSystem2D
typedef BaryCoordSystem< 3 > Go::BaryCoordSystem3D

Functions

template<typename T, int Dim>
Go::Array< T, Dim > Go::operator * (T d, const Go::Array< T, Dim > &v)
 The product of a vector and a scalar.
template<typename T, int Dim>
std::istream & Go::operator>> (std::istream &is, Go::Array< T, Dim > &v)
 Stream extraction for Array.
template<typename T, int Dim>
std::ostream & Go::operator<< (std::ostream &os, const Go::Array< T, Dim > &v)
 Stream insertion for Array.
template<class T, int Dim>
Array< T, Dim > Go::operator * (const Array< double, Dim > &a, const T b)
template<int Dim>
std::istream & std::operator>> (std::istream &is, Go::BaryCoordSystem< Dim > &bc)
 Read BaryCoordSystem from input stream.
template<int Dim>
std::ostream & std::operator<< (std::ostream &os, const Go::BaryCoordSystem< Dim > &bc)
 Write BaryCoordSystem to output stream.
double Go::binom (int n, int i)
 Computes the binomial coefficient: n! / (i! (n-i)!).
double Go::factorial (int n)
 computes n! (n factorial)
double Go::trinomial (int n, int i, int j)
 computes the trinomial coefficient: n! / (i! j! (n-i-j)!)
double Go::quadrinomial (int n, int i, int j, int k)
 computes the quadrinomial coefficient: n! / (i! j! k! (n-i-j-k)!)
std::istream & std::operator>> (std::istream &is, Go::BoundingBox &bbox)
 Read BoundingBox from input stream.
std::ostream & std::operator<< (std::ostream &os, const Go::BoundingBox &bbox)
 Write BoundingBox to output stream.
template<class Functor>
double Go::brent_minimize (const Functor &f, double a, double b, double c, double &parmin, const double rel_tolerance=std::sqrt(std::numeric_limits< double >::epsilon()))
double Go::curvatureRadius (const std::vector< Point > &der, std::vector< Point > &unitder)
 Given position, first and second derivative of a curve passing through a point, compute the unit tangent, curvature vector and curvature radius of this curve.
double Go::stepLenFromRadius (double radius, double aepsge)
 Computes the step length along a curve based on radius of curvature at a point on the curve, and an absolute tolerance.
double Go::tanLenFromRadius (double radius, double angle)
 To create the tangent length for interpolating a circular arc with an almost equi-oscillating Hermit qubic.
void Go::getHermiteData (const std::vector< Point > &der1, const std::vector< Point > &der2, double &parint, double &len1, double &len2)
 Given position, first and second derivative in both ends of an Hermite segment, compute parameter interval and tangent lengths in order to stay close to a circular segment.
template<class Functor>
void Go::minimise_conjugated_gradient (FunctionMinimizer< Functor > &dfmin)
 This is the algorithm for minimising a function taking multiple parameters, using the conjugated gradient method.
template<typename Functor>
void Go::trapezoidal (Functor &f, double a, double b, double &s, int n)
 This routine computes the n'th stage of refinement of an extended trapezoidal rule.
template<typename Functor>
double Go::simpsons_rule (Functor &f, double a, double b, const double eps=1.0e-6, const int max_iter=20)
 Routine to calculate the integral of the functor f from a to b using Simpson's rule.
template<typename Functor>
double Go::gaussian_quadrature (Functor &f, double a, double b)
 Routine to calculate the integral of the functor f from a to b using Gaussian quadrature with W=1 and N=10.
template<typename Functor2D>
double Go::simpsons_rule2D (Functor2D &f, double ax, double bx, double ay, double by, const double eps=1.0e-6, const int max_iter=20)
 Routine to integrate the two-dimensional functor f over the rectangle defined by ax, bx, ay and by.
template<typename Functor2D>
double Go::gaussian_quadrature2D (Functor2D &f, double ax, double bx, double ay, double by)
 Routine to integrate the two-dimensional functor f over the rectangle defined by ax, bx, ay and by.
template<typename SquareMatrix>
void Go::LUDecomp (SquareMatrix &mat, int num_rows, int *perm, bool &parity)
 LU decomposition algorithm, based on Crout's algorithm.
template<typename SquareMatrix, typename T>
void Go::LUsolveSystem (SquareMatrix &A, int num_unknowns, T *vec)
 Solve the system Ax = b for x, using LU decomposition of the matrix A.
template<typename SquareMatrix, typename T>
void Go::forwardSubstitution (const SquareMatrix &L, T *x, int num_unknowns)
 Using forward substitution to calculate x on the system Lx = b, where L is a lower triangular matrix with unitary diagonal.
template<typename SquareMatrix, typename T>
void Go::backwardSubstitution (const SquareMatrix &U, T *x, int num_unknowns)
 Using backward substitution to calculate x on the system Ux = b, where U is an upper triangular matrix with unitary diagonal.
template<typename SquareMatrix>
void Go::forwardSubstitution (const SquareMatrix &L, std::vector< double > *x, int num_unknowns)
 Using forward substitution to calculate x on the system Lx = b, where L is a lower triangular matrix with unitary diagonal.
template<typename SquareMatrix>
void Go::backwardSubstitution (const SquareMatrix &U, std::vector< double > *x, int num_unknowns)
 Using backward substitution to calculate x on the system Ux = b, where U is an upper triangular matrix with unitary diagonal.
template<typename T, int Dim>
std::ostream & Go::operator<< (std::ostream &os, const MatrixXD< T, Dim > &m)
 output operator
Point Go::operator * (double d, const Point &p)
 The product of a vector and a scalar.
std::istream & Go::operator>> (std::istream &is, Go::Point &v)
 Stream extraction for Point.
std::ostream & Go::operator<< (std::ostream &os, const Go::Point &v)
 Stream insertion for Point.
bool Go::operator< (const Point &p1, const Point &p2)
void Go::normalNoise (double *res, double mean_err, int num_samples)
 Gives a certain number of random samples drawn from the normal distribution.
void Go::uniformNoise (double *res, double lval, double uval, int num_samples)
 Gives a certain number of random samples drawn from the uniform distribution.
Rational Go::operator+ (const Rational &r1, const Rational r2)
Rational Go::operator- (const Rational &r1, const Rational r2)
Rational Go::operator * (const Rational &r1, const Rational r2)
Rational Go::operator/ (const Rational &r1, const Rational r2)
std::ostream & Go::operator<< (std::ostream &os, const Rational &p)
double Go::getCurrentTime ()
 Number of seconds since some (probably system-dependent) epoch.
void Go::systemSleep (double sleep_time)
 Sleep for sleep_time seconds.
template<typename T>
Go::determinantOf (const Array< T, 2 > *a)
 Calculates the determinant of a 2 x 2 matrix, represented in memory as a sequence of 2 Array s of length 2.
template<typename T>
Go::determinantOf (const Array< T, 3 > *a)
 Calculates the determinant of a 3 x 3 matrix, represented in memory as a sequence of 3 Array s of length 3.
template<typename T, int Dim>
Go::simplex_volume (const Array< T, Dim > *a)
 Computes the volume of a simplex consisting of (Dim+1) vertices embedded in Euclidean space of dimension (Dim).
template<typename T>
Go::area (const Array< T, 2 > *c)
 Computes the area of a 2-dimensional triangle.
template<typename T>
Go::area (const Array< T, 3 > *c)
 Computes the area of a 2-dimensional triangle.
template<typename T>
Go::volume (const Array< T, 3 > *c)
 Computes the volume of a 3D simplex (embedded i 3D space).
template<typename T>
Go::signed_area (const Array< T, 3 > *c, const Array< T, 3 > &normal)
 Computes the signed area of a triangle embedded in 3D space.
template<class T>
Array< T, Dim > Go::BaryCoordSystem::baryToCart (const Array< T, Dim+1 > &bary_pt) const
 Input is a barycentric point, output is the corresponding cartesian point.
template<class T>
Array< T, Dim+1 > Go::BaryCoordSystem::cartToBary (const Array< T, Dim > &cart_pt) const
 Input is a cartesian point, output is the corresponding barycentric point.
 Go::FunctionMinimizer::FunctionMinimizer (int num_param, const Functor &fun, const double *const seed, double tol=std::sqrt(std::numeric_limits< double >::epsilon()))
 Constructor.
double Go::FunctionMinimizer::minimize (const Point &dir, bool &hit_domain_edge)
 Move the current parameter point in order to minimize the function along an unoriented line.
double Go::FunctionMinimizer::linminBrent (const Point &dir, const double *bracket, const double *fval_brak)
 Move the current parameter point in order to minimize the function along a ray.
double Go::FunctionMinimizer::fval () const
 evaluating distance function at the current point of evaluation
double Go::FunctionMinimizer::fval (const Point &param) const
 evaluating the distance function for an arbitrary parameter value (not the parameter pair).
void Go::FunctionMinimizer::grad (Point &result) const
 evaluate the gradient (2 components) of the distance function at the currently set parameter pair of evaluation.
void Go::FunctionMinimizer::grad (const Point &param, Point &result) const
 evaluate the gradient (2 components) of the distance function for an arbitrary parameter value (not the current parameter pair).
void Go::FunctionMinimizer::moveUV (const Point &dir, double multiplier)
 movie the current parameter pair a certain distance ('multiplier') in a certain direction ('dir', which has 2 components)
void Go::MatrixXD::setToRotation (T angle, T x, T y, T z)
 Similar to glRotate(), sets the matrix to be a rotation of angle radians about the axis given by (x, y, z).
void Go::MatrixXD< double, 3 >::setToRotation (doubleangle, doublex, doubley, doublez)
 Similar to glRotate(), sets the matrix to be a rotation of angle radians about the axis given by (x, y, z).
void Go::MatrixXD::setToRotation (const Vector3D &p, const Vector3D &q)
 Sets the matrix to be a rotation that takes the point p to q.
void Go::MatrixXD< double, 3 >::setToRotation (const Vector3D &p, const Vector3D &q)
 Sets the matrix to be a rotation that takes the point p to q.

Function Documentation

template<typename T>
T Go::area ( const Array< T, 3 > *  c  )  [inline]

Computes the area of a 2-dimensional triangle.

Input is an array of corner points. Same function also exists for 2-dimensional triangles.

Definition at line 94 of file Volumes.h.

References Go::Array< T, Dim >::cross(), and Go::Array< T, Dim >::length().

template<typename T>
T Go::area ( const Array< T, 2 > *  c  )  [inline]

Computes the area of a 2-dimensional triangle.

Input is an array of corner points. Same function also exists for 3-dimensional triangles.

Definition at line 87 of file Volumes.h.

References Go::simplex_volume().

Referenced by Go::BaryCoordSystemTriangle3D::BaryCoordSystemTriangle3D().

template<typename SquareMatrix>
void Go::backwardSubstitution ( const SquareMatrix &  U,
std::vector< double > *  x,
int  num_unknowns 
)

Using backward substitution to calculate x on the system Ux = b, where U is an upper triangular matrix with unitary diagonal.

Parameters:
U The upper triangular matrix. The actually used class must support the operation [][] and return 'double'.
x At function invocation, x should point to a vector of doubles, containing the vector 'b'. On successful completion of the function, this vector will contain the solution for x.
num_unknowns The system size (number of unknowns).

Definition at line 188 of file LUDecomp_implementation.h.

template<typename SquareMatrix, typename T>
void Go::backwardSubstitution ( const SquareMatrix &  U,
T *  x,
int  num_unknowns 
)

Using backward substitution to calculate x on the system Ux = b, where U is an upper triangular matrix with unitary diagonal.

Parameters:
U The upper triangular matrix. The actually used class must support the operation [][] and return 'double'.
x At function invocation, x should point to an array of T, containing the vector 'b'. On successful completion of the function, this array will contain the solution for x.
num_unknowns The system size (number of unknowns).

Definition at line 174 of file LUDecomp_implementation.h.

Referenced by Go::LUsolveSystem().

template<typename T>
T Go::determinantOf ( const Array< T, 3 > *  a  )  [inline]

Calculates the determinant of a 3 x 3 matrix, represented in memory as a sequence of 3 Array s of length 3.

Same function also exists for 2 x 2 matrices.

Definition at line 61 of file Volumes.h.

template<typename T>
T Go::determinantOf ( const Array< T, 2 > *  a  )  [inline]

Calculates the determinant of a 2 x 2 matrix, represented in memory as a sequence of 2 Array s of length 2.

Same function also exists for 3 x 3 matrices.

Definition at line 51 of file Volumes.h.

Referenced by Go::simplex_volume().

template<typename SquareMatrix>
void Go::forwardSubstitution ( const SquareMatrix &  L,
std::vector< double > *  x,
int  num_unknowns 
)

Using forward substitution to calculate x on the system Lx = b, where L is a lower triangular matrix with unitary diagonal.

Parameters:
L The lower triangular matrix. The actually used class must support the operation [][] and return 'double'.
x At function invocation, x should point to a vector of doubles, containing the vector 'b'. On successful completion of the function, this vector will contain the solution for x.
num_unknowns The system size (number of unknowns).

Definition at line 159 of file LUDecomp_implementation.h.

template<typename SquareMatrix, typename T>
void Go::forwardSubstitution ( const SquareMatrix &  L,
T *  x,
int  num_unknowns 
)

Using forward substitution to calculate x on the system Lx = b, where L is a lower triangular matrix with unitary diagonal.

Parameters:
L The lower triangular matrix. The actually used class must support the operation [][] and return 'double'.
x At function invocation, x should point to an array of T, containing the vector 'b'. On successful completion of the function, this array will contain the solution for x.
num_unknowns The system size (number of unknowns).

Definition at line 147 of file LUDecomp_implementation.h.

Referenced by Go::LUsolveSystem().

template<class Functor>
Go::FunctionMinimizer< Functor >::FunctionMinimizer ( int  num_param,
const Functor &  fun,
const double *const   seed,
double  tol = std::sqrt(std::numeric_limits< double >::epsilon()) 
) [inline, inherited]

Constructor.

'minimization_tol' should in general be kept no lower than the square root of machine precision. 'seed' is a pointer to an array containing the start values for the search in the 'num_param' parameters.

Definition at line 164 of file GeneralFunctionMinimizer_implementation.h.

template<typename Functor2D>
double Go::gaussian_quadrature2D ( Functor2D &  f,
double  ax,
double  bx,
double  ay,
double  by 
)

Routine to integrate the two-dimensional functor f over the rectangle defined by ax, bx, ay and by.

Uses Gaussian quadrature.

Definition at line 168 of file Integration.h.

References Go::gaussian_quadrature().

template<class Functor>
double Go::FunctionMinimizer< Functor >::linminBrent ( const Point dir,
const double *  bracket,
const double *  fval_brak 
) [inline, inherited]

Move the current parameter point in order to minimize the function along a ray.

The ray is given by the current parameter point and the parameter 'dir'. In contrast to minimize(), linminBrent() requires that the minimum is bracketed. The method uses an algorithm inspired by Brent's method (with certain modifications). This method is rapid when the function is well-behaved. Otherwise it will use a slow but robust 'golden mean'-search. The current parameter point is moved to the minimum found along the given ray. Before calling this function, the minimum must be bracketed, and the brackets must be given to the function through the 'brackets' variable, and the corresponding function values through the 'fval_brak' variable. The function value at the new minimum is returned.

Definition at line 321 of file GeneralFunctionMinimizer_implementation.h.

References Go::FunctionMinimizer< Functor >::fval(), and Go::FunctionMinimizer< Functor >::moveUV().

Referenced by Go::brent_minimize(), and Go::FunctionMinimizer< Functor >::minimize().

template<typename SquareMatrix>
void Go::LUDecomp ( SquareMatrix &  mat,
int  num_rows,
int *  perm,
bool &  parity 
)

LU decomposition algorithm, based on Crout's algorithm.

Parameters:
mat The matrix to be decomposed. The actually used class must support the operation [][] and return 'double'.
num_rows Number of rows (which is also equal to the number of columns).
perm Should point to an n-sized array where the permutation is written.
parity Value upon function completion is 'true' if the number of row interchanges was pair. 'false' otherwise.

Definition at line 48 of file LUDecomp_implementation.h.

References Go::sum().

Referenced by Go::LUsolveSystem().

template<typename SquareMatrix, typename T>
void Go::LUsolveSystem ( SquareMatrix &  A,
int  num_unknowns,
T *  vec 
)

Solve the system Ax = b for x, using LU decomposition of the matrix A.

Upon successful completion of the function, the matrix A will be LU decomposed, and the solution x will be computed and stored at the memory area pointed to by the last argument.

Parameters:
A The system matrix A. The actually used class must support the operation [][] and return 'double'.
num_unknowns Number of unknowns in the equation system.
vec At function invocation, 'vec' should point to an array of T, containing the vector 'b'. On successful completion of the function, this array will contain the solution for x.

Definition at line 128 of file LUDecomp_implementation.h.

References Go::backwardSubstitution(), Go::forwardSubstitution(), and Go::LUDecomp().

template<class Functor>
void Go::minimise_conjugated_gradient ( FunctionMinimizer< Functor > &  dfmin  ) 

This is the algorithm for minimising a function taking multiple parameters, using the conjugated gradient method.

It is used in conjunction with the FunctionMinimizer class. Documentation for both can be found here: FunctionMinimizer

Definition at line 74 of file GeneralFunctionMinimizer_implementation.h.

References Go::FunctionMinimizer< Functor >::atMax(), Go::FunctionMinimizer< Functor >::atMin(), Go::FunctionMinimizer< Functor >::fval(), Go::FunctionMinimizer< Functor >::grad(), Go::Point::length2(), Go::FunctionMinimizer< Functor >::minimize(), and Go::FunctionMinimizer< Functor >::numPars().

template<class Functor>
double Go::FunctionMinimizer< Functor >::minimize ( const Point dir,
bool &  hit_domain_edge 
) [inline, inherited]

Move the current parameter point in order to minimize the function along an unoriented line.

The line is given by the current parameter point and the parameter 'dir'.

Parameters:
dir a direction parallel to the search line
Return values:
hit_domain_edge indicates whether the found minimum is on the boundary of the domain.

Definition at line 203 of file GeneralFunctionMinimizer_implementation.h.

References Go::FunctionMinimizer< Functor >::fval(), Go::FunctionMinimizer< Functor >::grad(), Go::FunctionMinimizer< Functor >::linminBrent(), Go::FunctionMinimizer< Functor >::moveUV(), and Go::FunctionMinimizer< Functor >::numPars().

Referenced by Go::minimise_conjugated_gradient().

void Go::normalNoise ( double *  res,
double  mean_err,
int  num_samples 
)

Gives a certain number of random samples drawn from the normal distribution.

Parameters:
res a pointer to the array where the resulting samples should be written.
mean_err the sigma parameter to the normal distribution.
num_samples the desired number of samples (should also be the size of the array pointed to by res.

template<typename T, int Dim>
std::ostream& Go::operator<< ( std::ostream &  os,
const MatrixXD< T, Dim > &  m 
) [inline]

output operator

Terminates the det() template recursion.

Definition at line 433 of file MatrixXD.h.

void Go::MatrixXD< double, 3 >::setToRotation ( const Vector3D p,
const Vector3D q 
) [inline, inherited]

Sets the matrix to be a rotation that takes the point p to q.

p and q are assumed to lie on the 2-sphere. Only makes sense for Dim == 3.

Definition at line 329 of file MatrixXD.h.

References Go::MatrixXD< T, Dim >::identity().

template<typename T, int Dim>
void Go::MatrixXD< T, Dim >::setToRotation ( const Vector3D p,
const Vector3D q 
) [inline, inherited]

Sets the matrix to be a rotation that takes the point p to q.

p and q are assumed to lie on the 2-sphere. Only makes sense for Dim == 3.

Definition at line 321 of file MatrixXD.h.

void Go::MatrixXD< double, 3 >::setToRotation ( double   angle,
double   x,
double   y,
double   z 
) [inline, inherited]

Similar to glRotate(), sets the matrix to be a rotation of angle radians about the axis given by (x, y, z).

(x, y, z) does not have to be normalized. Only makes sense for Dim == 3 or Dim == 4.

Definition at line 268 of file MatrixXD.h.

References Go::MatrixXD< T, Dim >::identity(), and Go::Array< T, Dim >::normalize().

template<typename T, int Dim>
void Go::MatrixXD< T, Dim >::setToRotation ( angle,
x,
y,
z 
) [inline, inherited]

Similar to glRotate(), sets the matrix to be a rotation of angle radians about the axis given by (x, y, z).

(x, y, z) does not have to be normalized. Only makes sense for Dim == 3 or Dim == 4.

Definition at line 261 of file MatrixXD.h.

Referenced by Go::MatrixXD< T, Dim >::setToRotation().

template<typename T>
T Go::signed_area ( const Array< T, 3 > *  c,
const Array< T, 3 > &  normal 
)

Computes the signed area of a triangle embedded in 3D space.

Input is an array of corner points and a normal to determine the sign.

Definition at line 113 of file Volumes.h.

References Go::Array< T, Dim >::cross(), and Go::Array< T, Dim >::length().

Referenced by Go::BaryCoordSystemTriangle3D::cartToBary().

template<typename Functor2D>
double Go::simpsons_rule2D ( Functor2D &  f,
double  ax,
double  bx,
double  ay,
double  by,
const double  eps = 1.0e-6,
const int  max_iter = 20 
)

Routine to integrate the two-dimensional functor f over the rectangle defined by ax, bx, ay and by.

Uses Simpson's rule.

Definition at line 156 of file Integration.h.

References Go::simpsons_rule().

template<typename Functor>
void Go::trapezoidal ( Functor &  f,
double  a,
double  b,
double &  s,
int  n 
)

This routine computes the n'th stage of refinement of an extended trapezoidal rule.

Used as a help routine for 'simpsons_rule', found further down in this file. NB: It only carries out the computation for ONE stage of the procedure, so if you want to use it to integrate a function, you must call it successively for n = 1, n = 2, ....

Definition at line 56 of file Integration.h.

References Go::sum().

Referenced by Go::simpsons_rule().

void Go::uniformNoise ( double *  res,
double  lval,
double  uval,
int  num_samples 
)

Gives a certain number of random samples drawn from the uniform distribution.

Parameters:
res a pointer to the array where the resulting samples should be written.
lval lower bound of the range from which the samples can take their values.
uval upper bound of the range from which the samples can take their values.
num_samples the desired number of samples (should also be the size of the array pointed to by res.


Generated on Mon Jun 11 14:48:18 2007 for GoTools Core Library by  doxygen 1.5.1