Basic Image AlgorithmS Library  2.8.0
 All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Files | Classes | Functions
Algorithms

Algorithms that solve Equation Systems or compute Singular values. More...

+ Collaboration diagram for Algorithms:

Files

file  Lapack.hh
 BIAS c++ wrapper class for lapack fortran77 functions.
 

Classes

class  BIAS::BSplineCurve
 this class is intended for drawing of quadratic and cubic regular B-Splines from given control points. More...
 
class  BIAS::FFT
 Wrapper to the fftw3 library. More...
 
class  BIAS::GaussHelmert
 see Foerstner, 2004, "Uncertainty And Projective Geometry" ! untested ! More...
 
class  BIAS::Histogram1D< T >
 Simple one dimensional histogram computation. More...
 
class  BIAS::Interpolator
 this class interpolates a function y=f(t) between given control points (the y-values) More...
 
class  BIAS::LaguerreSolver
 class encapsulating a laguerre solver for polynomials More...
 
class  BIAS::LeastSquaresBase
 Base class for linear least squares solvers. More...
 
class  BIAS::LeastSquaresLapack
 Linear least squares solver based on Lapack routines. More...
 
class  BIAS::LeastSquaresSVD
 Linear least squares solver based on singular value decomposition. More...
 
class  BIAS::Median1D< DataType >
 Computes the median and p-quantile of a vector. More...
 
class  BIAS::MonteCarloTransform
 monte carlo propagation of uncertainty More...
 
class  BIAS::MSAC< SolutionType >
 this class does something, but WHAT? Is it the M-Estimator SAmple Consensus ?? More...
 
class  BIAS::PolynomialSolve
 base class for solving polynomial equations More...
 
class  BIAS::PreemptiveRANSAC< SolutionType >
 Fast RANSAC after David Nister, "Preemptive RANSAC for Live Structure And Motion Estimation", Internation Conference on Computer Vision (ICCV) 2003. More...
 
class  BIAS::Random
 class for producing random numbers from different distributions More...
 
class  BIAS::RANSAC< SolutionType >
 Generic abstract base class for RANSAC (RANdom SAmple Consensus) estimators. More...
 
class  BIAS::RANSACPreKnowledge< SolutionType >
 Extension of RANSAC algorithm: user-supplied guesses are evaluated and refined using RANSAC criterions before actually starting the RANSAC algorithm. More...
 
class  BIAS::StandardDeviation
 Computes the standard deviation and mean of given values. More...
 
class  BIAS::SVD
 computes and holds the singular value decomposition of a rectangular (not necessarily quadratic) Matrix A More...
 
class  BIAS::SVD3x3
 singular value decomposition for 3x3 matrices More...
 
class  BIAS::UncertaintyTransformBase
 base class for all uncertainty transforms More...
 
class  BIAS::UnscentedTransform
 uses the unscented transformation to map a normal distribututed random variable using a nonlinear transformation More...
 
class  BIAS::VectorStatistics
 Provides methods to compute statistical values from vector<double>. More...
 

Functions

void BIAS::CheckJacobian (minpack_funcder_mn fun, void *p, const long int NumErrors, const Vector< double > &x, Vector< double > &error)
 Convenience function to check correctness of user supplied Jacobians J(x) needed by the analytical version of the Levenberg-Marquardt algorithm. More...
 
long int BIAS::ComputeCovariance (long int M, long int N, double *ResultingJacobian, int *ipvt, double &SumOfSquaredErrors, Matrix< double > &Cov)
 Compute covariance matrix from Levenberg-Marquardt resulting Jacobian matrix J(x) and permutation matrix P (see description for method LevenbergMarquardtExtended above or search in minpack/minpack-documentation.txt for IPVT). More...
 
long int BIAS::ComputeCovariance (const long int NumErrors, const long int NumParams, const Matrix< double > &Jac, const Vector< int > &Permutation, const double &SumOfSquaredErrors, Matrix< double > &Cov)
 Convenience wrapper for ComputeCovariance. More...
 
long int BIAS::ComputeJacobian (minpack_func_mn fun, void *p, const int NumErrors, const Vector< double > &x, Matrix< double > &Jac, const double gradientEpsilon=LM_DEF_EPSILON)
 Compute the Jacobian J(x) of function f at x. More...
 
BIASMathAlgo_EXPORT int Eigenproblem_quadratic_matrix (const BIAS::Matrix< double > &A, BIAS::Vector< double > &ret_EigenValuesReal, BIAS::Vector< double > &ret_EigenValuesImag, BIAS::Matrix< double > &eigenVectors)
 solve general eigenproblem for a general quadratix (n x n) matrix. More...
 
BIASMathAlgo_EXPORT int Eigenvalue_solve (const BIAS::Matrix< double > &A, BIAS::Vector< double > &wr, BIAS::Vector< double > &wi)
 solve unsymmetric eigenvalue problems (for a rectangular matrix) untested, grest More...
 
BIASMathAlgo_EXPORT
BIAS::Matrix< double > 
Fortran_Matrix_to_Matrix (const TNT::Fortran_Matrix< double > &FMat)
 converts the argument Fortran_Matrix to a Matrix and returns it. More...
 
BIASMathAlgo_EXPORT void Fortran_Matrix_to_Matrix (const TNT::Fortran_Matrix< double > &FMat, BIAS::Matrix< double > &res)
 converts the argument Fortran_Matrix to a Matrix and returns it. More...
 
BIASMathAlgo_EXPORT int General_singular_value_decomposition (char jobu, char jobvt, const BIAS::Matrix< double > &A, BIAS::Vector< double > &ret_S, BIAS::Matrix< double > &ret_U, BIAS::Matrix< double > &ret_VT)
 solve general eigenproblem. More...
 
BIASMathAlgo_EXPORT int General_singular_value_decomposition (const BIAS::Matrix< double > &A, BIAS::Vector< double > &ret_S, BIAS::Matrix< double > &ret_U, BIAS::Matrix< double > &ret_VT)
 
BIASMathAlgo_EXPORT int General_singular_value_decomposition (const BIAS::Matrix< double > &A, BIAS::Vector< double > &ret_S, BIAS::Matrix< double > &ret_VT)
 
BIASMathAlgo_EXPORT
BIAS::Matrix< double > 
generalised_eigenvalue_matrix_solve (BIAS::Matrix< double > &A_, BIAS::Matrix< double > &B_)
 
BIASMathAlgo_EXPORT int Lapack_Cholesky_SymmetricPositiveDefinit (const BIAS::Matrix< double > &A, BIAS::Matrix< double > &UL, const char ul='U')
 Coputes the Cholesky decomposition of a symmetric positive definit matrix A. More...
 
BIASMathAlgo_EXPORT
BIAS::Vector< double > 
Lapack_LLS_QR_linear_solve (const BIAS::Matrix< double > &A, const BIAS::Vector< double > &b, int &res)
 linear least squares solves |Ax-b|=min res=0 on success, anything else rudimentary tested (woelk) More...
 
BIASMathAlgo_EXPORT
BIAS::Vector< double > 
Lapack_LLS_QR_linear_solve (const BIAS::Matrix< double > &A, const BIAS::Vector< double > &b)
 as above but ignores success info (res) More...
 
BIASMathAlgo_EXPORT
BIAS::Vector< double > 
Lapack_LU_linear_solve (const BIAS::Matrix< double > &A, const BIAS::Vector< double > &b)
 solve linear equations using LU factorization More...
 
BIASMathAlgo_EXPORT
BIAS::Vector< double > 
Lapack_WLLS_solve (const BIAS::Matrix< double > &A, const BIAS::Vector< double > &b, const BIAS::Matrix< double > &B, int &res)
 weighted linear least squares solves |B^-1(b-Ax)|=min res=0 on success, anything else rudimentary tested More...
 
long int BIAS::LevenbergMarquardt (minpack_funcder_mn fun, void *p, long int m, long int n, Vector< double > &InitialGuess, Vector< double > &Result, double Tolerance=MINPACK_DEF_TOLERANCE)
 Solves an overdetermined system f(x) = 0 of m non-linear equations with n parameters using the Levenberg-Marquardt algorithm. More...
 
long int BIAS::LevenbergMarquardt (minpack_func_mn fun, void *p, long int M, long int N, Vector< double > &InitialGuess, Vector< double > &Result, double Tolerance=MINPACK_DEF_TOLERANCE, long int MaxIter=LM_DEF_MAX_ITER)
 Solves an overdetermined system f(x) = 0 of m non-linear equations with n parameters using the Levenberg-Marquardt algorithm. More...
 
long int BIAS::LevenbergMarquardtExtended (minpack_funcder_mn fun, void *p, long int m, long int n, Vector< double > &InitialGuess, Vector< double > &Result, double Tolerance=MINPACK_DEF_TOLERANCE, long int MaxIter=LM_DEF_MAX_ITER)
 Solves an overdetermined system f(x) = 0 of m non-linear equations with n parameters using the extended Levenberg-Marquardt algorithm. More...
 
long int BIAS::LevenbergMarquardtExtended (minpack_func_mn fun, void *p, long int M, long int N, Vector< double > &InitialGuess, Vector< double > &Result, double EpsilonFun=LM_DEF_EPSILON, double Tolerance=MINPACK_DEF_TOLERANCE, long int MaxIter=LM_DEF_MAX_ITER, double *ResultingJacobian=NULL, int *ipvt=NULL, double *SumOfSquaredErrors=NULL)
 Solves an overdetermined system f(x) = 0 of m non-linear equations with n parameters using the Levenberg-Marquardt algorithm. More...
 
BIASMathAlgo_EXPORT double MahalanobisDistance (const BIAS::Matrix< double > &Sigma, const BIAS::Vector< double > &d)
 computes squared Mahalanobis distance V^T Sigma^-1 V efficiently using cholesky decomposition and exploitation of symmetry More...
 
BIASMathAlgo_EXPORT int Packed_symmetric_eigenvalue_solve (long int N, const BIAS::Vector< double > &A, BIAS::Vector< double > &eigVals, std::vector< BIAS::Vector< double > > &eigVecs, const bool upperTriangle=true)
 Solve symmetric eigenvalue problem for matrix in packed storage. More...
 
long int BIAS::Powell (minpack_func_nn fun, void *p, Vector< double > &InitialGuess, Vector< double > &Result, double Tolerance=MINPACK_DEF_TOLERANCE)
 Solve a system f(x) = 0 of n non-linear equations with n unknowns using the Powell hybrid method. More...
 
long int BIAS::PowellExtended (minpack_func_nn fun, void *p, Vector< double > &InitialGuess, Vector< double > &Result, double EpsilonFun=LM_DEF_EPSILON, double Tolerance=MINPACK_DEF_TOLERANCE)
 Solve a system f(n) = 0 of n non-linear equations with n unknowns using the Powell hybrid method. More...
 
void BIAS::QRFrac (const Matrix< double > &A, Matrix< double > &Q, Matrix< double > &R)
 Compute the QR decomposition of input matrix A using Householder transformations without column pivoting. More...
 
BIASMathAlgo_EXPORT double SquaredMahalanobisDistance (const BIAS::Matrix< double > &Sigma, const BIAS::Vector< double > &d)
 computes Mahalanobis distance V^T Sigma^-1 V efficiently using cholesky decomposition and exploitation of symmetry More...
 
BIASMathAlgo_EXPORT
BIAS::Vector< double > 
Upper_symmetric_eigenvalue_solve (const BIAS::Matrix< double > &A)
 Solve symmetric eigenvalue problem (eigenvalues only) More...
 
BIASMathAlgo_EXPORT int Upper_symmetric_eigenvalue_solve (const BIAS::Matrix< double > &A, BIAS::Vector< double > &eigenvalues, std::vector< BIAS::Vector< double > > &eigenvectors)
 Compute eigenvalues and eigenvectors of a symmetric real square matrix. More...
 

Detailed Description

Algorithms that solve Equation Systems or compute Singular values.

C++ Wrappers for important Lapack fortran77 functions.

culultaive beta distribution

C++ Wrappers for gnu scientific library (GSL) p(x) dx = {(a+b) (a) (b)} x^{a-1} (1-x)^{b-1} dx for 0 <= x <= 1.

C++ Wrappers for gnu scientific library (GSL) 1/Gamma(x)

C++ Wrappers for gnu scientific library (GSL) (x) = ^ dt t^{x-1} (-t) The gamma function, x cannot be greater than 171.0 and x cannot be a negative integer.

Computes the inverse culmulative probability at x of a chi square distribution with deg_freedom degrees of freedom C++ Wrappers for gnu scientific library (GSL)

Computes the culmulative probability at x of a chi square distribution with deg_freedom degrees of freedom C++ Wrappers for gnu scientific library (GSL)

Computes the probability density at x of a chi square distribution with deg_freedom degrees of freedom C++ Wrappers for gnu scientific library (GSL)

Wrapper to Lapack and Minpack Functions.

Author
woelk 07/2008 (c) www.vision-n.de
woelk 02/2006

For numerical reasons, a+b must be less than GSL_SF_GAMMA_XMAX (171.0)

Author
woelk 02/2006
woelk 09/2006

For numerical reasons, a+b must be less equal GSL_SF_GAMMA_XMAX (171.0), however when a==b the gamma distribution is approximated by a gaussian distribution. When a equals b, the beta distribution is symmetric about x=0.5. It is approximated by a gaussian distribution for very large a and b

Author
woelk 02/2006

Functions from Lapack, e.g. solving systems of linear equations or computing eigenvalues For further information have a look here http://www.netlib.org/lapack/lug/.

Function Documentation

BIASMathAlgo_EXPORT void BIAS::CheckJacobian ( minpack_funcder_mn  fun,
void *  p,
const long int  NumErrors,
const Vector< double > &  x,
Vector< double > &  error 
)

Convenience function to check correctness of user supplied Jacobians J(x) needed by the analytical version of the Levenberg-Marquardt algorithm.

First, the user defined Jacobian J(x) of f at x is calculated by calling fun with iflag = 2. Then the minpack function chkder is used to evaluate the Jacobian.

Warning
CheckJacobian does not perform reliably if cancellation or rounding errors cause a severe loss of significance in the evaluation of a function!
Parameters
funThe target function fun is supposed to have the same interface as described above in LevenbergMarquardt (version WITH Jacobian computation):

int fun(void *p, int m, int n, const double *x, double *fvec, double *fjac, int ldfjac, int iflag)

Parameters
ppoints to an object with additional data used for evaluating J(x)
NumErrorsis the number m of residuals of f(x) (i.e. the number of rows of J(x))
xis a vector of size n that defines where the Jacobian J(x) is evaluated at
errorThe function returns an error measure for each row of the Jacobian J(x) in the output vector error of size NumErrors.

On output, error contains measures of correctness of the respective gradients. If there is no severe loss of significance, then if error[i] is 1.0, the i-th gradient is correct, while if error[i] is 0.0, the i-th gradient is incorrect. For values of error[i] between 0.0 and 1.0, the categorization is less certain. In general, a value of error[i] > 0.5 indicates that the i-th gradient is probably correct, while a value of error[i] < 0.5 indicates that the i-th gradient is probably incorrect.

Author
woelk 07/2006

Definition at line 406 of file Minpack.cpp.

References BIAS::Matrix< T >::GetData(), BIAS::Vector< T >::GetData(), TNT::Vector< T >::newsize(), and TNT::Vector< T >::size().

BIASMathAlgo_EXPORT long int BIAS::ComputeCovariance ( long int  M,
long int  N,
double *  ResultingJacobian,
int *  ipvt,
double &  SumOfSquaredErrors,
Matrix< double > &  Cov 
)

Compute covariance matrix from Levenberg-Marquardt resulting Jacobian matrix J(x) and permutation matrix P (see description for method LevenbergMarquardtExtended above or search in minpack/minpack-documentation.txt for IPVT).

Parameters
MNumber of residuals (= rows of Jacobian)
NNumber of parameters (= columns of Jacobian)
ResultingJacobianPointer to array of size M*N containing the resulting Jacobian J(x_res) computed from call to ExtendedLevenbergMarquardt.
ipvtPointer to integer array of size N describing the permutation matrix P computed from call to ExtendedLevenbergMarquardt (i.e. J(x_res)*P = Q*R).
SumOfSquaredErrorsSum of squared errors computed by ExtendedLevenbergMarquardt
CovOutput matrix to store covariance in
Author
woelk 06/2006

Definition at line 289 of file Minpack.cpp.

References BIAS::Matrix< T >::GetData(), and BIAS::Vector< T >::GetData().

BIASMathAlgo_EXPORT long int BIAS::ComputeCovariance ( const long int  NumErrors,
const long int  NumParams,
const Matrix< double > &  Jac,
const Vector< int > &  Permutation,
const double &  SumOfSquaredErrors,
Matrix< double > &  Cov 
)

Convenience wrapper for ComputeCovariance.

Compute covariance matrix from Levenberg-Marquardt resulting Jacobian matrix J(x) and permutation matrix P (see description for method LevenbergMarquardtExtended above or search in minpack/minpack-documentation.txt for IPVT).

Parameters
NumErrorsNumber of residuals (= rows of Jacobian)
NumParamsNumber of parameters (= columns of Jacobian)
JacMatrix of size NumErrors x NumParams containing the resulting Jacobian J(x_res) computed from ExtendedLevenbergMarquardt.
PermutationInteger vector of size NumParams describing the permutation matrix P computed from ExtendedLevenbergMarquardt (i.e. J(x_res)*P = Q*R).
SumOfSquaredErrorsSum of squared errors computed by ExtendedLevenbergMarquardt
CovOutput matrix to store covariance int
Author
koeser

Definition at line 304 of file Minpack.cpp.

References BIAS::Equal(), BIAS::Matrix< T >::GetData(), BIAS::Matrix< T >::Mult(), BIAS::Matrix< T >::MultiplyIP(), TNT::Matrix< T >::newsize(), BIAS::Matrix< T >::SetZero(), and TMP.

BIASMathAlgo_EXPORT long int BIAS::ComputeJacobian ( minpack_func_mn  fun,
void *  p,
const int  NumErrors,
const Vector< double > &  x,
Matrix< double > &  Jac,
const double  gradientEpsilon = LM_DEF_EPSILON 
)

Compute the Jacobian J(x) of function f at x.

Parameters
funspecifies the target function f(x) to compute the Jacobian J(x) for
ppoints to an object with additional data used for evaluating J(x)
xis a vector of size n that defines where the Jacobian J(x) is evaluated at
Jacis an output matrix to which the Jacobian J(x) of function f evaluated at x is written to
Note
Wrapper for fdjac2
See Also
minpack/fdjac2.c
Author
woelk 07/2007

Definition at line 443 of file Minpack.cpp.

References BIAS::Matrix< T >::GetData(), BIAS::Vector< T >::GetData(), TNT::Vector< T >::size(), and BIAS::Matrix< T >::Transpose().

Referenced by BIAS::LevenbergMarquardtBase::LM_ComputeWithoutJacobian().

BIASMathAlgo_EXPORT int Eigenproblem_quadratic_matrix ( const BIAS::Matrix< double > &  A,
BIAS::Vector< double > &  ret_EigenValuesReal,
BIAS::Vector< double > &  ret_EigenValuesImag,
BIAS::Matrix< double > &  eigenVectors 
)

solve general eigenproblem for a general quadratix (n x n) matrix.

computes eigenvalues and eigenvectors of A. (Jan Woetzel), 03/2002

Definition at line 626 of file Lapack.cpp.

References TNT::Fortran_Matrix< T >::begin(), BIAS::Matrix< T >::GetData(), TNT::Matrix< T >::newsize(), TNT::Fortran_Matrix< T >::num_cols(), TNT::Matrix< T >::num_cols(), TNT::Fortran_Matrix< T >::num_rows(), TNT::Matrix< T >::num_rows(), TNT::Matrix< T >::size(), and TNT::transpose().

BIASMathAlgo_EXPORT int Eigenvalue_solve ( const BIAS::Matrix< double > &  A,
BIAS::Vector< double > &  wr,
BIAS::Vector< double > &  wi 
)

solve unsymmetric eigenvalue problems (for a rectangular matrix) untested, grest

Definition at line 580 of file Lapack.cpp.

References TNT::Matrix< T >::num_cols(), TNT::Fortran_Matrix< T >::num_rows(), TNT::Matrix< T >::num_rows(), and TNT::transpose().

BIASMathAlgo_EXPORT BIAS::Matrix<double> Fortran_Matrix_to_Matrix ( const TNT::Fortran_Matrix< double > &  FMat)

converts the argument Fortran_Matrix to a Matrix and returns it.

it mainly tranposes the data array because Fortran-Matrices different row/col-order

Author
Jan Woetzel bug fixed by Daniel Grest (04/2003)

Definition at line 169 of file Lapack.cpp.

References TNT::Fortran_Matrix< T >::num_cols(), and TNT::Fortran_Matrix< T >::num_rows().

Referenced by General_singular_value_decomposition().

BIASMathAlgo_EXPORT void Fortran_Matrix_to_Matrix ( const TNT::Fortran_Matrix< double > &  FMat,
BIAS::Matrix< double > &  res 
)

converts the argument Fortran_Matrix to a Matrix and returns it.

it mainly tranposes the data array because Fortran-Matrices different row/col-order, same as above but now copying necessary

Author
Daniel Grest beta April 2003

Definition at line 190 of file Lapack.cpp.

References BIAS::Matrix< T >::GetData(), TNT::Matrix< T >::newsize(), TNT::Fortran_Matrix< T >::num_cols(), TNT::Matrix< T >::num_cols(), TNT::Fortran_Matrix< T >::num_rows(), and TNT::Matrix< T >::num_rows().

BIASMathAlgo_EXPORT int General_singular_value_decomposition ( char  jobu,
char  jobvt,
const BIAS::Matrix< double > &  A,
BIAS::Vector< double > &  ret_S,
BIAS::Matrix< double > &  ret_U,
BIAS::Matrix< double > &  ret_VT 
)

solve general eigenproblem.

computes singular value decomposition of a rectangular m x n matrix A calls extern liblapack routine. better use SVD.hh ret_S vector with singular values of A with s(i) >= s(i+1) ret_U m x m orthogonal/unitary matrix U ret_VT n x n orthogonal/unitary matrx transpose(V) Jan Woetzel 03/2002, changed for Matrix by Daniel Grest, July 2002 woelk: added jobu and jobvt jobu: 'A' - calculate U 'N' - do not calculate U (much faster) not tested see man dgesvd

Definition at line 688 of file Lapack.cpp.

References Fortran_Matrix_to_Matrix(), TNT::Vector< T >::newsize(), TNT::Matrix< T >::num_cols(), TNT::Matrix< T >::num_rows(), and TNT::transpose().

Referenced by BIAS::HMatrixEstimation::Compute(), BIAS::HMatrixEstimation::ComputeAffine(), BIAS::SVD::General_Eigenproblem_GeneralMatrix_Lapack(), General_singular_value_decomposition(), BIAS::LeastSquaresSVD::Solve(), and BIAS::LeastSquaresSVD::WeightedSolve().

BIASMathAlgo_EXPORT int General_singular_value_decomposition ( const BIAS::Matrix< double > &  A,
BIAS::Vector< double > &  ret_S,
BIAS::Matrix< double > &  ret_U,
BIAS::Matrix< double > &  ret_VT 
)
inline

Definition at line 221 of file Lapack.hh.

References General_singular_value_decomposition().

BIASMathAlgo_EXPORT int General_singular_value_decomposition ( const BIAS::Matrix< double > &  A,
BIAS::Vector< double > &  ret_S,
BIAS::Matrix< double > &  ret_VT 
)
inline

Definition at line 230 of file Lapack.hh.

References General_singular_value_decomposition().

BIASMathAlgo_EXPORT BIAS::Matrix<double> generalised_eigenvalue_matrix_solve ( BIAS::Matrix< double > &  A_,
BIAS::Matrix< double > &  B_ 
)
BIASMathAlgo_EXPORT int Lapack_Cholesky_SymmetricPositiveDefinit ( const BIAS::Matrix< double > &  A,
BIAS::Matrix< double > &  UL,
const char  ul = 'U' 
)

Coputes the Cholesky decomposition of a symmetric positive definit matrix A.

When ul='U', the factorization has the form A = UL^T * UL where UL is a upper triangular matrix. When ul='L',the factorization has the form A = UL * UL^T where UL is a lower triangular matrix. Wrapper for lapacks dpotf2.

Author
woelk 09/2007

Definition at line 212 of file Lapack.cpp.

References BIAS::Matrix< T >::GetCols(), BIAS::Matrix< T >::GetRows(), TNT::Matrix< T >::newsize(), TNT::Matrix< T >::num_cols(), TNT::Fortran_Matrix< T >::num_rows(), TNT::Matrix< T >::num_rows(), BIAS::Matrix< T >::SetZero(), and TNT::transpose().

Referenced by BIAS::UnscentedTransform::ComputeSigmaPoints_(), and SquaredMahalanobisDistance().

BIASMathAlgo_EXPORT BIAS::Vector<double> Lapack_LLS_QR_linear_solve ( const BIAS::Matrix< double > &  A,
const BIAS::Vector< double > &  b,
int &  res 
)
BIASMathAlgo_EXPORT BIAS::Vector<double> Lapack_LLS_QR_linear_solve ( const BIAS::Matrix< double > &  A,
const BIAS::Vector< double > &  b 
)
inline

as above but ignores success info (res)

Definition at line 73 of file Lapack.hh.

References Lapack_LLS_QR_linear_solve().

BIASMathAlgo_EXPORT BIAS::Vector<double> Lapack_LU_linear_solve ( const BIAS::Matrix< double > &  A,
const BIAS::Vector< double > &  b 
)
BIASMathAlgo_EXPORT BIAS::Vector<double> Lapack_WLLS_solve ( const BIAS::Matrix< double > &  A,
const BIAS::Vector< double > &  b,
const BIAS::Matrix< double > &  B,
int &  res 
)

weighted linear least squares solves |B^-1(b-Ax)|=min res=0 on success, anything else rudimentary tested

Author
woelk 04 2003

Definition at line 376 of file Lapack.cpp.

References BIAS::Vector< T >::GetData(), TNT::Matrix< T >::num_cols(), TNT::Matrix< T >::num_rows(), TNT::Vector< T >::size(), and TNT::transpose().

Referenced by BIAS::LeastSquaresLapack::WeightedSolve().

BIASMathAlgo_EXPORT long int BIAS::LevenbergMarquardt ( minpack_funcder_mn  fun,
void *  p,
long int  m,
long int  n,
Vector< double > &  InitialGuess,
Vector< double > &  Result,
double  Tolerance = MINPACK_DEF_TOLERANCE 
)

Solves an overdetermined system f(x) = 0 of m non-linear equations with n parameters using the Levenberg-Marquardt algorithm.

A least squares solution is used. Jacobian J(x) of target function f(x) is given analytically.

Parameters
funDefines the target function f(x) (resp. the Jacobian matrix J(x)) which is supposed to have the following interface:

int fun(void *p, int m, int n, const double *x, double *fvec, double *fjac, int ldfjac, int iflag)

  • p points to an object with additional data used for evaluating f(x)
  • m is number of residuals (size of array fvec)
  • n is number of parameters (size of array x)
  • x points to an array of size n containing the input parameter values
  • fvec points to an array of size m where the output residuals f(x) are written to
  • fjac points to an array of size m*n where the elements of the Jacobian matrix J(x) are written to
  • ldfjac specifies the leading dimension used to describe the Jacobian matrix J(x) by the array fjac (if ldfjac = m then J(x) is written column-wise into array fjac).
  • iflag is a flag that specifies if the function should compute the residuals f(x) and write them to fvec (when iflag == 1), or if the Jacobian matrix J(x) should be computed and written to fjac (when iflag == 2). The value of iflag should not be changed!
Parameters
mNumber of residuals f(x) with m>=n
nNumber of variables x with m>=n
InitialGuessStarting point x_0 for optimization (input, size must be n)
ResultSolution x_res after optimization (output, size must be n)
ToleranceNon-negative threshold used to determine convergence.
Returns
0 determines success in optimization, other return values are:
  • return = -1: Improper input parameters
  • return = 0: A termination condition became true.
  • return = 5: Number of calls to fun has reached or exceeded maxfev iterations.
  • return = 6: Tolerance ftol is too small. No further reduction in sum of residual squares possible.
  • return = 7: Tolerance xtol is too small. No further improvement in approximate solution x is possible.
  • return = 8: Tolerance gtol is too small. Vector fvec is orthogonal to the columns of the Jacobian to machine precision.
See Also
minpack/lmder1.f
Author
woelk

Definition at line 97 of file Minpack.cpp.

References BIAS::Vector< T >::GetData().

Referenced by BIAS::PMatrixEstimation::AutoCalib_(), BIAS::PolynomialSolve::FitPolynomial(), BIAS::ProjectionParametersSphericalFast::FitPolynomial(), BIAS::ProjectionParametersSpherical::FitPolynomialWithoutLinearMonomial_(), BIAS::FMatrixEstimation::GoldStandard(), BIAS::HMatrixEstimation::Optimize(), and BIAS::FMatrixEstimation::Optimize().

BIASMathAlgo_EXPORT long int BIAS::LevenbergMarquardt ( minpack_func_mn  fun,
void *  p,
long int  M,
long int  N,
Vector< double > &  InitialGuess,
Vector< double > &  Result,
double  Tolerance = MINPACK_DEF_TOLERANCE,
long int  MaxIter = LM_DEF_MAX_ITER 
)

Solves an overdetermined system f(x) = 0 of m non-linear equations with n parameters using the Levenberg-Marquardt algorithm.

A least squares solution is used. The Jacobian J(x) of target function f(x) is calculated numerically by forward-difference approximation and must not be provided analytically.

Parameters
funDefines the target function f(x) which is supposed to have the following interface:

int fun(void *p, int m, int n, const double *x, double *fvec, int iflag)

  • p points to an object with additional data used for evaluating f(x)
  • m is number of residuals (size of array fvec)
  • n is number of parameters (size of array x)
  • x points to an array of size n containing the input parameter values
  • fvec points to an array of size m where the output residuals f(x) are written to
  • iflag can be ignores here
Parameters
mNumber of residuals f(x) with m>=n
nNumber of variables x with m>=n
InitialGuessStarting point x_0 for optimization (input, size must be n)
ResultSolution x_res after optimization (output, size must be n)
ToleranceNon-negative threshold used to determine convergence.
MaxIterThis parameter is ignored and should be removed!
Returns
0 determines success in optimization, for other return see description for method LevenbergMarquardt above.
See Also
minpack/lmdif1.f
Author
woelk

Definition at line 184 of file Minpack.cpp.

References BIAS::Vector< T >::GetData().

BIASMathAlgo_EXPORT long int BIAS::LevenbergMarquardtExtended ( minpack_funcder_mn  fun,
void *  p,
long int  m,
long int  n,
Vector< double > &  InitialGuess,
Vector< double > &  Result,
double  Tolerance = MINPACK_DEF_TOLERANCE,
long int  MaxIter = LM_DEF_MAX_ITER 
)

Solves an overdetermined system f(x) = 0 of m non-linear equations with n parameters using the extended Levenberg-Marquardt algorithm.

A least squares solution is used. Jacobian J(x) of target function f(x) is given analytically.

Parameters
funDefines the target function f(x) (resp. the Jacobian matrix J(x)) which is supposed to have the following interface:

int fun(void *p, int m, int n, const double *x, double *fvec, double *fjac, int ldfjac, int iflag)

  • p points to an object with additional data used for evaluating f(x)
  • m is number of residuals (size of array fvec)
  • n is number of parameters (size of array x)
  • x points to an array of size n containing the input parameter values
  • fvec points to an array of size m where the output residuals f(x) are written to
  • fjac points to an array of size m*n where the elements of the Jacobian matrix J(x) are written to
  • ldfjac specifies the leading dimension used to describe the Jacobian matrix J(x) by the array fjac (if ldfjac = m then J(x) is written column-wise into array fjac).
  • iflag is a flag that specifies if the function should compute the residuals f(x) and write them to fvec (when iflag == 1), or if the Jacobian matrix J(x) should be computed and written to fjac (when iflag == 2). The value of iflag should not be changed!
Parameters
mNumber of residuals f(x) with m>=n
nNumber of variables x with m>=n
InitialGuessStarting point x_0 for optimization (input, size must be n)
ResultSolution x_res after optimization (output, size must be n)
ToleranceNon-negative threshold used to determine convergence.
MaxIterMax. number of evaluations of target function f(x). Note that the complete numerical approximation of the Jacobian J(x) is always done as a whole before this value is checked. Set MaxIter to -1 if you want to use the default number of max. iterations from the fortran code (= 200*(n+1)).
Returns
0 determines success in optimization, for other return see description for method LevenbergMarquardt above.
See Also
minpack/lmder.f
Author
woelk

Definition at line 131 of file Minpack.cpp.

References BIAS::Vector< T >::GetData().

Referenced by BIAS::RMatrix::GetAverageRMatrix(), BIAS::LevenbergMarquardtBase::LM_Compute(), and BIAS::LevenbergMarquardtBase::LM_ComputeWithoutJacobian().

BIASMathAlgo_EXPORT long int BIAS::LevenbergMarquardtExtended ( minpack_func_mn  fun,
void *  p,
long int  M,
long int  N,
Vector< double > &  InitialGuess,
Vector< double > &  Result,
double  EpsilonFun = LM_DEF_EPSILON,
double  Tolerance = MINPACK_DEF_TOLERANCE,
long int  MaxIter = LM_DEF_MAX_ITER,
double *  ResultingJacobian = NULL,
int *  ipvt = NULL,
double *  SumOfSquaredErrors = NULL 
)

Solves an overdetermined system f(x) = 0 of m non-linear equations with n parameters using the Levenberg-Marquardt algorithm.

A least squares solution is used. The Jacobian J(x) of target function f(x) is calculated numerically by forward-difference approximation and must not be provided analytically.

Parameters
funDefines the target function f(x) which is supposed to have the following interface:

int fun(void *p, int m, int n, const double *x, double *fvec, int iflag)

  • p points to an object with additional data used for evaluating f(x)
  • m is number of residuals (size of array fvec)
  • n is number of parameters (size of array x)
  • x points to an array of size n containing the input parameter values
  • fvec points to an array of size m where the output residuals f(x) are written to
  • iflag can be ignores here
Parameters
mNumber of residuals f(x) with m>=n
nNumber of variables x with m>=n
InitialGuessStarting point x_0 for optimization (input, size must be n)
ResultSolution x_res after optimization (output, size must be n)
EpsilonFunNon-negative step width used to approximate the Jacobian numerically.
ToleranceNon-negative threshold used to determine convergence.
MaxIterMax. number of evaluations of target function f(x) (-1: default).
ResultingJacobianPointer to array of size m*n to write the Jacobian J(x_res) evaluated numerically at the solution x_res. Pass NULL if it should be omitted.
ipvtIf Jacobian J(x_res) should be returned in ResultingJacobian, you must specify here a pointer to an integer array of size n to write the permutation matrix P to. Otherwise pass NULL.

P is computed so that J(x_res)*P = Q*R, where J(x_res) is the final calculated Jacobian, Q is an orthogonal matrix, and R is an upper triangular matrix with diagonal elements of nonincreasing magnitude. Column j of P is column ipvt[j] of the identity matrix. Search minpack/minpack-documentation.txt for IPVT for a detailed description of this parameter!

Returns
0 determines success in optimization, for other return see description for method LevenbergMarquardt above.
See Also
minpack/lmdif.f
Author
woelk

Definition at line 218 of file Minpack.cpp.

References BIAS::Vector< T >::GetData().

BIASMathAlgo_EXPORT double MahalanobisDistance ( const BIAS::Matrix< double > &  Sigma,
const BIAS::Vector< double > &  d 
)

computes squared Mahalanobis distance V^T Sigma^-1 V efficiently using cholesky decomposition and exploitation of symmetry

   Sigma must be square and symmetric positive definite, same dim as V

   About 10 times faster (for 6x6) than V^T svd.invert(Sigma) V,
   but it still has potential for optimization, see doc in code.
Author
koeser 11/2007

Definition at line 780 of file Lapack.cpp.

References SquaredMahalanobisDistance().

BIASMathAlgo_EXPORT int Packed_symmetric_eigenvalue_solve ( long int  N,
const BIAS::Vector< double > &  A,
BIAS::Vector< double > &  eigVals,
std::vector< BIAS::Vector< double > > &  eigVecs,
const bool  upperTriangle = true 
)

Solve symmetric eigenvalue problem for matrix in packed storage.

Parameters
[in]NNumber of rows/cols of N x N matrix A
[in]AData vector of A in packed storage of length N*(N+1)/2. If upper triangle matrix is given, vector contains matrix rows above the diagonal (i.e. A11, A12, A13, ..., A1N, A22, A23, ...), otherweise vector contains matrix columns below the diagonal (i.e. A11, A21, A31, ..., AN1, A22, A32, ...).
[out]eigValsReturns eigenvalues in ascending order
[out]eigVecsReturn eigenvectors for eigenvalues
[in]upperTriangleSpecifies if packed vector for matrix A is given as upper (row-wise) or lower (column-wise) triangle.
Note
Wrapper for LAPACK routine dspev.
Author
esquivel 08/2013

Definition at line 535 of file Lapack.cpp.

References TNT::Vector< T >::begin(), TNT::Vector< T >::newsize(), BIAS::Vector< T >::SetZero(), and BIAS::Vector< T >::Size().

BIASMathAlgo_EXPORT long int BIAS::Powell ( minpack_func_nn  fun,
void *  p,
Vector< double > &  InitialGuess,
Vector< double > &  Result,
double  Tolerance = MINPACK_DEF_TOLERANCE 
)

Solve a system f(x) = 0 of n non-linear equations with n unknowns using the Powell hybrid method.

The Jacobian matrix J(x) is calculated numerically by forward-difference approximation. The error is assumed to be within the range of machine precision.

Parameters
funis the name of the user-supplied subroutine which calculates the target function f at x and returns the values in array fvec. fun must be declared e.g. in the user calling program and have the following interface:

int fun(void *p, int n, const double *x, double *fvec, int iflag)

  • p points to an object with additional data used for evaluating f(x)
  • n is a positive integer input variable set to the number of residuals and parameters (i.e. size of array x and fvec).
  • x is the parameter array of length n. On input, x must contain an initial estimate of the solution vector. On output, x contains the final estimate of the solution vector x_res.
  • fvec is an output array of length n which contains the function values f(x) evaluated at the input parameter vector x.
  • iflag should not be changed by fun unless the user wants to terminate execution, in this case set iflag to a negative integer.
Parameters
Toleranceis a non-negative input variable. Termination occurs when the algorithm estimates that the relative error between x and the solution x_res is at most Tolerance.
Returns
0 if everything worked fine. If the user has terminated execution by setting iflag to < 0 in fun, the value of iflag is returned. Otherwise, return is:
  • return = -1: Improper input parameters
  • return = 0: Algorithm estimated that the relative error between x and the solution is within tolerance range.
  • return = 2: Number of calls to fun has reached or exceeded 200*(n+1).
  • return = 3: Tolerance is too small. No further improvement in the approximate solution x is possible.
  • return = 4: Iteration is not making good progress.
See Also
minpack/hybrd1.f
Author
woelk
Examples:
ExamplePowell.cpp.

Definition at line 35 of file Minpack.cpp.

References BIAS::Vector< T >::GetData(), and BIAS::Vector< T >::Size().

Referenced by BIAS::PMatrixEstimation::ComputeFromFQuasiEuklid(), and BIAS::PolynomialSolve::NonLinearRefine().

BIASMathAlgo_EXPORT long int BIAS::PowellExtended ( minpack_func_nn  fun,
void *  p,
Vector< double > &  InitialGuess,
Vector< double > &  Result,
double  EpsilonFun = LM_DEF_EPSILON,
double  Tolerance = MINPACK_DEF_TOLERANCE 
)

Solve a system f(n) = 0 of n non-linear equations with n unknowns using the Powell hybrid method.

The Jacobian matrix J(x) is approximated numerically by forward-difference approximation.

See description of method Powell above for further parameter documentation.

Parameters
EpsilonFunNon-negative step width used to approximate the Jacobian numerically.
See Also
minpack/hybrd.f
Author
woelk

Definition at line 56 of file Minpack.cpp.

References BIAS::Vector< T >::GetData(), and BIAS::Vector< T >::Size().

BIASMathAlgo_EXPORT void BIAS::QRFrac ( const Matrix< double > &  A,
Matrix< double > &  Q,
Matrix< double > &  R 
)

Compute the QR decomposition of input matrix A using Householder transformations without column pivoting.

Parameters
[in]Ais the input matrix which will be QR decomposed.
[out]Qis the resulting Q matrix which is orthogonal (A = Q*R)
[out]Ris the resulting R matrix which is upper trapezoidal (A = Q*R)
Note
Wrapper for qrfac
See Also
minpack/qrfac.c

Definition at line 473 of file Minpack.cpp.

References BIAS::Matrix< T >::GetCols(), BIAS::Matrix< T >::GetRows(), TNT::Matrix< T >::newsize(), BIAS::Vector< T >::OuterProduct(), BIAS::Matrix< T >::SetIdentity(), BIAS::Vector< T >::SetZero(), BIAS::Matrix< T >::SetZero(), and BIAS::Matrix< T >::Transpose().

BIASMathAlgo_EXPORT double SquaredMahalanobisDistance ( const BIAS::Matrix< double > &  Sigma,
const BIAS::Vector< double > &  d 
)

computes Mahalanobis distance V^T Sigma^-1 V efficiently using cholesky decomposition and exploitation of symmetry

   Sigma must be square and symmetric positive definite, same dim as V

   About 10 times faster (for 6x6) than V^T svd.invert(Sigma) V,
   but it still has potential for optimization, see doc in code.
Author
koeser 11/2007

Definition at line 787 of file Lapack.cpp.

References Lapack_Cholesky_SymmetricPositiveDefinit(), Lapack_LU_linear_solve(), BIAS::Matrix< T >::MakeSymmetric(), BIAS::Matrix< T >::NormFrobenius(), TNT::Matrix< T >::num_rows(), and BIAS::Vector< T >::Size().

Referenced by MahalanobisDistance().

BIASMathAlgo_EXPORT BIAS::Vector<double> Upper_symmetric_eigenvalue_solve ( const BIAS::Matrix< double > &  A)

Solve symmetric eigenvalue problem (eigenvalues only)

Note
Wrapper for LAPACK routine dsyev. untested
Author
grest

Definition at line 438 of file Lapack.cpp.

References TNT::Vector< T >::begin(), TNT::Fortran_Matrix< T >::begin(), TNT::Matrix< T >::num_cols(), TNT::Fortran_Matrix< T >::num_rows(), TNT::Matrix< T >::num_rows(), and TNT::transpose().

Referenced by BIAS::MonteCarloTransform::GetSamples_(), BIAS::SVD::Sqrt(), and BIAS::SVD::SqrtT().

BIASMathAlgo_EXPORT int Upper_symmetric_eigenvalue_solve ( const BIAS::Matrix< double > &  A,
BIAS::Vector< double > &  eigenvalues,
std::vector< BIAS::Vector< double > > &  eigenvectors 
)

Compute eigenvalues and eigenvectors of a symmetric real square matrix.

Note
Wrapper for LAPACK routine dsyevd.
Author
woelk 11/2006

Definition at line 477 of file Lapack.cpp.

References TNT::Vector< T >::begin(), TNT::Vector< T >::newsize(), TNT::Matrix< T >::num_cols(), TNT::Fortran_Matrix< T >::num_rows(), TNT::Matrix< T >::num_rows(), and TNT::transpose().