Woolz Image Processing  Version 1.7.5
AlgMatrix

Files

file  AlgMatrix.c
 Matrix allocation and maintenance functions.
 
file  AlgMatrixCG.c
 Conjugate Gradient iterative method with preconditioning for the solution of linear systems with the form \(\mathbf{A} \mathbf{x} = \mathbf{b}\). A must be a symmetric postive definite matrix, i.e. \({\mathbf{x}}^T \mathbf{A} \mathbf{x} < 0\), \(\forall \mathbf{x} \not= \mathbf{0}\), \(\mathbf{x} \in R^n\).
 
file  AlgMatrixGauss.c
 Provides a function for solving matrix equations of the form: A.x = b for x using Gaussian elimination with partial pivoting.
 
file  AlgMatrixLSQR.c
 Provides functions for solving matrix equations using LSQR. This software is based on lsqr.c, a C version of LSQR derived by James W. Howse jhows.nosp@m.e@la.nosp@m.nl.go.nosp@m.v from the Fortran 77 implementation of C. C. Paige and M. A. Saunders. In most cases the extensive comments from Howse's lsqr.c have been preserved with little change.
 
file  AlgMatrixLU.c
 Provides functions for solving matrix equations of the form: A.x = b for x, inverting a matrix, calculating the determinant of a matrix and performing LU decomposition.
 
file  AlgMatrixMath.c
 Functions for basic arithmatic with matricies.
 
file  AlgMatrixRSEigen.c
 Functions to find the eigenvalues and eigenvectors of a real symmetric matrix.
 
file  AlgMatrixRSTDiag.c
 Reduces a real symmetric matrix to symmetric tridiagonal form by orthogonal similarity transformation and construction of the right operator of the reduction.
 
file  AlgMatrixSV.c
 Provides functions for singular value decomposition.
 
file  AlgMatrixTDiagQLI.c
 Determines the eigenvalues and eigenvectors of a real symmetric tridiagonal matrix using implicit shifts.
 

Functions

AlgMatrix AlgMatrixNew (AlgMatrixType aType, size_t nR, size_t nC, size_t nE, double tol, AlgError *dstErr)
 Allocates a new matrix or the requested type. More...
 
AlgMatrixRectAlgMatrixRectNew (size_t nR, size_t nC, AlgError *dstErr)
 Allocates a new rectangular matrix. More...
 
AlgMatrixSymAlgMatrixSymNew (size_t nN, AlgError *dstErr)
 Allocates a new symmetric matrix. More...
 
AlgMatrixLLRAlgMatrixLLRNew (size_t nR, size_t nC, size_t nE, double tol, AlgError *dstErr)
 Allocates a new linked list row matrix. More...
 
AlgMatrixLLREAlgMatrixLLRENew (AlgMatrixLLR *mat)
 Gets a linked list row matrix entry from the matrix. The return value may be NULL if none are available and the entries need to be expanded sing AlgMatrixLLRExpand(). More...
 
void AlgMatrixLLREFree (AlgMatrixLLR *mat, AlgMatrixLLRE *p)
 Returns a linked list row matrix entry to the free stack of the matrix. More...
 
AlgError AlgMatrixWriteAscii (AlgMatrix mat, FILE *fP)
 Writes a matrix in numeric ASCI format to the given file file. The rows are on separate lines and the columns of each row are white space seperated. More...
 
AlgMatrix AlgMatrixReadAscii (AlgMatrixType mType, double tol, FILE *fP, const char *fSep, size_t recMax, AlgError *dstErr)
 Reads a matrix of the requested type from an ASCII file. If the type is ALG_MATRIX_SYM then only the first half of each row will be used, although all must be given. More...
 
AlgError AlgMatrixSet (AlgMatrix mat, size_t row, size_t col, double val)
 
AlgError AlgMatrixLLRCopyInPlace (AlgMatrixLLR *aM, AlgMatrixLLR *bM)
 Copy the second LLR matrix to the first. More...
 
AlgError AlgMatrixLLRSet (AlgMatrixLLR *mat, size_t row, size_t col, double val)
 This can have one of three actions: More...
 
double AlgMatrixValue (AlgMatrix mat, size_t row, size_t col)
 Returns the value in the matrix at the given coordinates. More...
 
double AlgMatrixLLRValue (AlgMatrixLLR *mat, size_t row, size_t col)
 Returns the value in the matrix at the given coordinates. More...
 
AlgError AlgMatrixLLRExpand (AlgMatrixLLR *mat, size_t nE)
 Ensures that there are at least the requested number of free linked list row matrix entries available. More...
 
void AlgMatrixZero (AlgMatrix mat)
 Sets all elements of the matrix to zero. More...
 
void AlgMatrixRectZero (AlgMatrixRect *mat)
 Sets all elements of the matrix to zero. More...
 
void AlgMatrixSymZero (AlgMatrixSym *mat)
 Sets all elements of the matrix to zero. More...
 
void AlgMatrixLLRZero (AlgMatrixLLR *mat)
 Sets all elements of the matrix to zero. More...
 
AlgError AlgMatrixSetAll (AlgMatrix mat, double val)
 Sets all elements of the given matrix to the given value. More...
 
void AlgMatrixRectSetAll (AlgMatrixRect *mat, double val)
 Sets all elements of the given rectangular matrix to the given value. More...
 
void AlgMatrixSymSetAll (AlgMatrixSym *mat, double val)
 Sets all elements of the given symmetric matrix to the given value. More...
 
AlgError AlgMatrixLLRSetAll (AlgMatrixLLR *mat, double val)
 Sets all elements of the given linked list row matrix to the given value. This is not an efficient use of a linked list row matrix since all values will be allocate and set. More...
 
void AlgMatrixLLRERemove (AlgMatrixLLR *mat, size_t row, size_t col)
 Removes the entry at the given coordinates. More...
 
AlgError AlgMatrixCGSolve (AlgMatrix aM, double *xV, double *bV, AlgMatrix wM, void(*pFn)(void *, AlgMatrix, double *, double *), void *pDat, double tol, int maxItr, double *dstTol, int *dstItr)
 Conjugate Gradient iterative method with preconditioning for the solution of linear systems with the form \(\mathbf{A} \mathbf{x} = \mathbf{b}\). \(\mathbf{A}\) must be a symmetric postive definite matrix, i.e. \({\mathbf{x}}^T \mathbf{A} \mathbf{x} < 0\), \(\forall \mathbf{x} \not= \mathbf{0}\), \(\mathbf{x} \in R^n\). Convergence is tested using: \( \frac{\| \mathbf{b} - mathbf{A} \mathbf{x} \|} {\|\mathbf{b}\|} < \delta\). If the preconditioning function pFn is non NULL then it is called, passing the preconditioning data pDat, as: (*pFn)(void *pDat, double **aM, double *r, double *z, int n) to solve \(\mathbf{A} \mathbf{z} = \mathbf{r}\) for \(\mathbf{z}\) with the solution overwriting the initial contents of z. More...
 
AlgError AlgMatrixGaussSolve (AlgMatrix aMat, double *xMat)
 Solves the matrix equation A.x = b for x by Gaussian elimination with partial pivoting. Matrix A is the matrix of coefficients. More...
 
AlgError AlgMatrixSolveLSQR (AlgMatrix aM, double *bV, double *xV, double damping, double relErrA, double relErrB, long maxItr, long condLim, int *dstTerm, long *dstItr, double *dstFNorm, double *dstCondN, double *dstResNorm, double *dstResNormA, double *dstNormX)
 This function finds a solution x to the following problems: More...
 
AlgError AlgMatrixLUSolveRaw3 (double **aM, double *bV, int bSz)
 Solves the matrix equation A.x = b for x, where A is a 3x3 libAlc double array matrix. On return the matrix A is overwritten with its LU decomposition and the column vector b is overwritten with the solution column matrix x. More...
 
AlgError AlgMatrixLUSolveRaw4 (double **aM, double *bV, int bSz)
 Solves the matrix equation A.x = b for x, where A is a 4x4 libAlc double array matrix. On return the matrix A is overwritten with its LU decomposition and the column vector b is overwritten with the solution column matrix x. More...
 
AlgError AlgMatrixLUSolve (AlgMatrix aM, double *bV, int bSz)
 Solves the matrix equation A.x = b for x, where A is a square matrix. On return the matrix A is overwritten with its LU decomposition and the column vector b is overwritten with the solution column matrix x. More...
 
AlgError AlgMatrixLUSolveRaw (double **aM, int aSz, double *bV, int bSz)
 Solves the matrix equation A.x = b for x, where A is a square libAlc double array matrix. On return the matrix A is overwritten with its LU decomposition and the column vector b is overwritten with the solution column matrix x. More...
 
AlgError AlgMatrixLUInvertRaw3 (double **aM)
 Calculates the inverse of a 3x3 libAlc double array matrix. More...
 
AlgError AlgMatrixLUInvertRaw4 (double **aM)
 Calculates the inverse of a 4x4 libAlc double array matrix. More...
 
AlgError AlgMatrixLUInvert (AlgMatrix aM)
 Calculates the inverse of a square matrix. More...
 
AlgError AlgMatrixLUInvertRaw (double **aM, int aSz)
 Calculates the inverse of a square matrix. More...
 
AlgError AlgMatrixLUDetermRaw3 (double **aM, double *det)
 Calculates the determinant of a 3x3 double libAlc array matrix. The matrix is overwitten with its LU decomposition on return. More...
 
AlgError AlgMatrixLUDetermRaw4 (double **aM, double *det)
 Calculates the determinant of a 4x4 double libAlc array matrix. The matrix is overwitten with its LU decomposition on return. More...
 
AlgError AlgMatrixLUDeterm (AlgMatrix aM, double *det)
 Calculates the determinant of a square matrix. The matrix is overwitten with its LU decomposition on return. More...
 
AlgError AlgMatrixLUDetermRaw (double **aM, int aSz, double *det)
 Calculates the determinant of a square double libAlc array matrix. The matrix is overwitten with its LU decomposition on return. More...
 
AlgError AlgMatrixLUDecomp (AlgMatrix aM, int *iV, double *evenOdd)
 Replaces the given matrix with the LU decomposition of a row-wise permutation of itself. The given index vector is used to record the permutation effected by the partial pivoting. More...
 
AlgError AlgMatrixLUDecompRaw (double **aM, int aSz, int *iV, double *evenOdd)
 Replaces the given matrix with the LU decomposition of a row-wise permutation of itself. The given index vector is used to record the permutation effected by the partial pivoting. More...
 
AlgError AlgMatrixLUBackSub (AlgMatrix aM, int *iV, double *bV)
 Solves the set of of linear equations A.x = b where A is input as its LU decomposition determined with AlgMatrixLUDecomp() of a row-wise permutation of itself. The given index vector is used to record the permutation effected by the partial pivoting. More...
 
AlgError AlgMatrixLUBackSubRaw (double **aM, int aSz, int *iV, double *bV)
 Solves the set of of linear equations A.x = b where A is input as its LU decomposition determined with AlgMatrixLUDecompRaw() of a row-wise permutation of itself. The given index vector is used to record the permutation effected by the partial pivoting. More...
 
void AlgMatrixAdd (AlgMatrix aM, AlgMatrix bM, AlgMatrix cM)
 Computes the sum of two matricies and returns the result in a third supplied matrix:

\[ \mathbf{A} = \mathbf{B} + \mathbf{C} \]

The dimensions of the all three matricies must be nR, nC. It is safe to supply the same matrix as any combination of aM, bM and cM. More...

 
void AlgMatrixSub (AlgMatrix aM, AlgMatrix bM, AlgMatrix cM)
 Subtracts on matrix from another and returns the result in a third supplied matrix:

\[ \mathbf{A} = \mathbf{B} - \mathbf{C} \]

The dimensions of the all three matricies must be nR, nC. It is safe to supply the same matrix as any combination of aM, bM and cM. More...

 
void AlgMatrixMul (AlgMatrix aM, AlgMatrix bM, AlgMatrix cM)
 Computes the product of two matricies and returns the result in a third supplied matrix:

\[ \mathbf{A} = \mathbf{B} \mathbf{C} \]

The dimensions of the result matrix (aM) must be cR, bC. All the matrices must be valid and of the same type except in the case of multiplying tow symmetric matrices when the result matrix must be rectangular (because in general the product of two symmetric matrices is not symmetric). More...

 
double AlgMatrixTrace (AlgMatrix aM)
 Computes the trace of the given matrix. More...
 
void AlgMatrixTranspose (AlgMatrix aM, AlgMatrix bM)
 Computes the transpose of the given matrix:

\[ \mathbf{A} = \mathbf{B^T} \]

aM = transpose(bM). The dimensions of the result matrix (aM) must be aR == bC, aC == bR. More...

 
void AlgMatrixCopy (AlgMatrix aM, AlgMatrix bM)
 Copies the values of the matrix bM to the result matric aM:

\[ \mathbf{A} = \mathbf{B} \]

. More...

 
void AlgMatrixScale (AlgMatrix aM, AlgMatrix bM, double sv)
 Multiplies the given matrix by the given scalar:

\[ \mathbf{A} = s \mathbf{B} \]

. More...

 
void AlgMatrixScaleAdd (AlgMatrix aM, AlgMatrix bM, AlgMatrix cM, double sv)
 Multiplies the a matrix by a scalar and then adds another matrix:

\[ \mathbf{A} = \mathbf{B} + s \mathbf{C}. \]

All the matrices must have the same dimensions. More...

 
void AlgMatrixScalar (AlgMatrix aM, double sv)
 Sets the elements of the given square matrix so that it is a scalar matrix:

\[ \mathbf{A} = s \mathbf{I} \]

. More...

 
void AlgMatrixVectorMul (double *aV, AlgMatrix bM, double *cV)
 Multiplies the matrix \(\mathbf{B}\) by the vector \(\mathbf{c}\):

\[ \mathbf{a} = \mathbf{B} \mathbf{c} \]

. More...

 
void AlgMatrixVectorMulAdd (double *aV, AlgMatrix bM, double *cV, double *dV)
 Multiplies the matrix \(\mathbf{B}\) by the vector \(\mathbf{c}\) and adds the vector \(\mathbf{d}\):

\[ \mathbf{a} = \mathbf{B} \mathbf{c} + \mathbf{d} \]

. More...

 
void AlgMatrixVectorMulWAdd (double *aV, AlgMatrix bM, double *cV, double *dV, double s, double t)
 Multiplies the matrix \(\mathbf{B}\) by the vector \(\mathbf{c}\) and adds the vector \(\mathbf{d}\) using the given weights \(s\) and \(t\):

\[ \mathbf{a} = s \mathbf{B} \mathbf{c} + t \mathbf{d} \]

. More...

 
void AlgMatrixTVectorMul (double *aV, AlgMatrix bM, double *cV)
 Multiplies the transpose of matrix \(\mathbf{B}\) by the vector \(\mathbf{c}\):

\[ \mathbf{a} = \mathbf{B}^T \mathbf{c} \]

. More...

 
void AlgMatrixTVectorMulAdd (double *aV, AlgMatrix bM, double *cV, double *dV)
 Multiplies the transpose of matrix \(\mathbf{B}\) by the vector \(\mathbf{c}\):

\[ \mathbf{a} = \mathbf{B}^T \mathbf{c} + \mathbf{d} \]

. More...

 
AlgError AlgMatrixRawInv2x2 (double *a00, double *a01, double *a10, double *a11)
 Inverts a raw 2 x 2 matrix. All matrix values are passed using pointers and their values are set to those of the inverse matrix on return. If the given matrix is singular the ALG_ERR_MATRIX_SINGULAR error code is returned and no values are changed. More...
 
AlgError AlgMatrixRawInv3x3 (double *a00, double *a01, double *a02, double *a10, double *a11, double *a12, double *a20, double *a21, double *a22)
 Inverts a raw 2 x 2 matrix. All matrix values are passed using pointers and their values are set to those of the inverse matrix on return. If the given matrix is singular the ALG_ERR_MATRIX_SINGULAR error code is returned and no values are changed. More...
 
AlgError AlgMatrixRSEigen (AlgMatrix aM, double *vM, int reqEV)
 Determines the eigenvalues and eigenvectors of a real symmetric matrix by calling AlgMatrixRSTDiag() to create a tridiagonal symmetric matrix and then AlgMatrixTDiagQLI() to compute its eigenvalues and eigenvectors. The eigenvectors and eigenvalues are returned in descending eigenvalue order. For efficiency, the eigenvectors should only be computed if required. More...
 
AlgError AlgMatrixRSTDiag (AlgMatrix aMat, double *dVec, double *oVec)
 An implementation of Householder's alorithm which reduces a real aSz x aSz symmetric matrix to symmetric tridiagonal form by aSz - 2 orthogonal similarity transformations. This is a translation of the functions/subroutines 'tred2' in the fortran EISPACK library and the Numerical Recipies book, both of these are descended from the algol procedure tred2 in Wilkinson J.H and Reinsch C. Linear Algebra, Volume II Handbook for Automatic Computation. Springer-Verlag, 1971. See Numerical Recipies in C: The Art of Scientific Computing. Cambridge University Press, 1992. More...
 
AlgError AlgMatrixSVSolve (AlgMatrix aMat, double *bVec, double tol, int *dstIC)
 Solves the matrix equation A.x = b for x, where A is a matrix with at least as many rows as columns. On return the matrix A is overwritten by the matrix U in the singular value decomposition: A = U.W.V'. More...
 
AlgError AlgMatrixSVDecomp (AlgMatrix aMat, double *wVec, AlgMatrix vMat)
 Performs singular value decomposition of the given matrix ( \(A\) ) and computes two additional matricies ( \(U\) ) and ( \(V\) ) such that:

\[ A = U W V^T \]

where \(V^T\) is the transpose of \(V\). The code for AlgMatrixSVDecomp was derived from: Numerical Recipies function svdcmp(), EISPACK subroutine SVD(), CERN subroutine SVD() and ACM algorithm 358. The singular values ( \(W\) ) and singular vectors ( \(V\) ) are returned unsorted. See AlgMatrixSVSolve() for a usage example. More...

 
AlgError AlgMatrixSVBackSub (AlgMatrix uMat, double *wVec, AlgMatrix vMat, double *bVec)
 Solves the set of of linear equations A.x = b where A is input as its singular value decomposition in the three matricies U, W and V, as returned by AlgMatrixSVDecomp(). The code for AlgMatrixSVBackSub was derived from: Numerical Recipies function svbksb(). More...
 
AlgError AlgMatrixTDiagQLI (double *dM, double *oM, int aSz, AlgMatrix zM)
 Determines the eigenvalues and eigenvectors of a real symmetric tridiagonal matrix using the QL algorithm with implicit shifts. This is a based on the function 'tqli' in Numerical Recipies in C: The Art of Scientific Computiung. Cambridge University Press, 1992. More...
 

Detailed Description

Function Documentation

AlgMatrix AlgMatrixNew ( AlgMatrixType  aType,
size_t  nR,
size_t  nC,
size_t  nE,
double  tol,
AlgError dstErr 
)

Allocates a new matrix or the requested type.

Returns
New matrix or matrix with core NULL on error.
Parameters
aTypeMatrix type.
nRNumber of rows.
nCNumber of columns.
nENumber of list entries to allocate, only used for linked list row matrices, may be zero.
tolMatrix tollerance value, only used for linked list row matrices.
dstErrDestination error pointer, may be NULL.

References ALG_ERR_MATRIX_TYPE, ALG_ERR_NONE, ALG_MATRIX_LLR, ALG_MATRIX_RECT, ALG_MATRIX_SYM, AlgMatrixLLRNew(), AlgMatrixRectNew(), AlgMatrixSymNew(), _AlgMatrix::core, _AlgMatrix::llr, _AlgMatrix::rect, and _AlgMatrix::sym.

Referenced by WlzKrigReallocBuffers2D(), and WlzKrigReallocBuffers3D().

AlgMatrixSym* AlgMatrixSymNew ( size_t  nN,
AlgError dstErr 
)

Allocates a new symmetric matrix.

Returns
New symmetric matrix or NULL on error.
Parameters
nNNumber of rows and columns.
dstErrDestination error pointer, may be NULL.

References AlcCalloc(), AlcFree(), AlcMalloc(), ALG_ERR_MALLOC, ALG_ERR_NONE, ALG_MATRIX_SYM, _AlgMatrixSym::array, _AlgMatrixSym::maxN, _AlgMatrixSym::nC, _AlgMatrixSym::nR, and _AlgMatrixSym::type.

Referenced by AlgMatrixNew(), and AlgMatrixReadAscii().

AlgMatrixLLR* AlgMatrixLLRNew ( size_t  nR,
size_t  nC,
size_t  nE,
double  tol,
AlgError dstErr 
)

Allocates a new linked list row matrix.

Returns
New linked list row matrix or NULL on error.
Parameters
nRNumber of rows.
nCNumber of columns.
nENumber of list entries to allocate, may be zero.
tolMatrix tollerance value.
dstErrDestination error pointer, may be NULL.

References AlcCalloc(), ALG_ERR_MALLOC, ALG_ERR_NONE, ALG_MATRIX_LLR, AlgMatrixLLRExpand(), AlgMatrixLLRFree(), _AlgMatrixLLR::nC, _AlgMatrixLLR::nR, _AlgMatrixLLR::tbl, and _AlgMatrixLLR::type.

Referenced by AlgMatrixNew(), AlgMatrixReadAscii(), and WlzCMeshCompSurfMap().

AlgMatrixLLRE* AlgMatrixLLRENew ( AlgMatrixLLR mat)

Gets a linked list row matrix entry from the matrix. The return value may be NULL if none are available and the entries need to be expanded sing AlgMatrixLLRExpand().

Returns
Linked list row matrix entry or NULL if none available.
Note
Call this function rather than manipulating the free entries directly.
Parameters
matLinked list row matrix.

References _AlgMatrixLLR::freeStk, _AlgMatrixLLR::numEnt, and _AlgMatrixLLRE::nxt.

Referenced by AlgMatrixLLRCopyInPlace(), AlgMatrixLLRSet(), AlgMatrixReadAscii(), and AlgMatrixScalar().

void AlgMatrixLLREFree ( AlgMatrixLLR mat,
AlgMatrixLLRE p 
)

Returns a linked list row matrix entry to the free stack of the matrix.

Note
Call this function rather than manipulating the free entries directly.
Parameters
matLinked list row matrix.
pLinked list row matrix entry.

References _AlgMatrixLLR::freeStk, _AlgMatrixLLR::numEnt, and _AlgMatrixLLRE::nxt.

Referenced by AlgMatrixLLRERemove(), and AlgMatrixLLRExpand().

AlgError AlgMatrixWriteAscii ( AlgMatrix  mat,
FILE *  fP 
)

Writes a matrix in numeric ASCI format to the given file file. The rows are on separate lines and the columns of each row are white space seperated.

Returns
Alg error code.
Parameters
matGiven matrix.
fPOutput file pointer.

References ALG_ERR_MATRIX_TYPE, ALG_ERR_NONE, ALG_ERR_WRITE, ALG_MATRIX_LLR, ALG_MATRIX_RECT, ALG_MATRIX_SYM, _AlgMatrixRect::array, _AlgMatrixSym::array, _AlgMatrixLLRE::col, _AlgMatrix::core, _AlgMatrix::llr, _AlgMatrixRect::nC, _AlgMatrixSym::nC, _AlgMatrixLLR::nC, _AlgMatrixRect::nR, _AlgMatrixSym::nR, _AlgMatrixLLR::nR, _AlgMatrixLLRE::nxt, _AlgMatrix::rect, _AlgMatrix::sym, _AlgMatrixLLR::tbl, _AlgMatrixCore::type, and _AlgMatrixLLRE::val.

Referenced by WlzBasisFnScalarMOS3DFromCPts(), and WlzCMeshCompSurfMap().

AlgMatrix AlgMatrixReadAscii ( AlgMatrixType  mType,
double  tol,
FILE *  fP,
const char *  fSep,
size_t  recMax,
AlgError dstErr 
)

Reads a matrix of the requested type from an ASCII file. If the type is ALG_MATRIX_SYM then only the first half of each row will be used, although all must be given.

Returns
New matrix read from the file. Matrix core member will be NULL on error.
Parameters
mTypeRequired type of matrix.
tolTolerance value for ALG_MATRIX_LLR matrices (unused for other types).
fPInput file pointer.
fSepField separator string containing possible field separator characters.
recMaxMaximum number of characters per input record.
dstErrDestination error pointer, may be NULL.

References ALC_ER_ALLOC, ALC_ER_NONE, ALC_ER_READ, AlcCalloc(), AlcVecReadDouble2Asci(), AlcVectorFree(), AlcVectorItemGet(), AlcVectorSetArray1D(), AlcVectorToArray2D(), ALG_ERR_FUNC, ALG_ERR_MALLOC, ALG_ERR_MATRIX_TYPE, ALG_ERR_NONE, ALG_ERR_READ, ALG_MATRIX_LLR, ALG_MATRIX_RECT, ALG_MATRIX_SYM, AlgMatrixFree(), AlgMatrixLLRENew(), AlgMatrixLLRExpand(), AlgMatrixLLRNew(), AlgMatrixSymNew(), _AlgMatrixRect::array, _AlgMatrixSym::array, _AlgMatrixLLRE::col, _AlgMatrix::core, _AlgMatrix::llr, _AlgMatrixRect::nC, _AlgMatrixRect::nR, _AlgMatrixLLRE::nxt, _AlgMatrix::rect, _AlgMatrix::sym, _AlgMatrixLLR::tbl, _AlgMatrixRect::type, and _AlgMatrixLLRE::val.

AlgError AlgMatrixSet ( AlgMatrix  mat,
size_t  row,
size_t  col,
double  val 
)
Returns
Alg error code.
            Errors can occur because of a memory allocation
            failure in ALG_MATRIX_LLR matrices or an invalid
            matrix type.
Parameters
matGiven matrix.
rowRow coordinate.
colColumn coordinate.
valMatrix value.

References ALG_ERR_NONE, ALG_MATRIX_LLR, ALG_MATRIX_RECT, ALG_MATRIX_SYM, AlgMatrixLLRSet(), _AlgMatrixRect::array, _AlgMatrixSym::array, _AlgMatrix::core, _AlgMatrix::llr, _AlgMatrix::rect, _AlgMatrix::sym, and _AlgMatrixCore::type.

Referenced by WlzCMeshCompSurfMap().

AlgError AlgMatrixLLRCopyInPlace ( AlgMatrixLLR aM,
AlgMatrixLLR bM 
)

Copy the second LLR matrix to the first.

Returns
Alg error code.
Parameters
aMFirst given matrix.
bMSecond given matrix.

References ALG_ERR_NONE, AlgMatrixLLRENew(), AlgMatrixLLRExpand(), AlgMatrixLLRZero(), _AlgMatrixLLRE::col, _AlgMatrixLLR::nR, _AlgMatrixLLR::numEnt, _AlgMatrixLLRE::nxt, _AlgMatrixLLR::tbl, and _AlgMatrixLLRE::val.

Referenced by AlgMatrixCopy().

AlgError AlgMatrixLLRSet ( AlgMatrixLLR mat,
size_t  row,
size_t  col,
double  val 
)

This can have one of three actions:

Returns
Alg error code.
            1. If no value exists at the given coordinates then
            the new value is added with a new entry.
            2. If a value exists at the given coordinate and the
            given value is non-zero then the entry value is
            replaced with the given value.
            3. If a value exists at the given coordinate and the
            given value is zero then the corresponding entry is
            removed.
            Errors can only occur because of a memory allocation
            failure. These can be avoided by preallocating
            sufficient entries using AlgMatrixLLRExpand().
Parameters
matLinked list row matrix.
rowRow coordinate.
colColumn coordinate.
valMatrix value.

References ALG_ERR_NONE, AlgMatrixLLRENew(), AlgMatrixLLRERemove(), AlgMatrixLLRExpand(), _AlgMatrixLLRE::col, _AlgMatrixLLR::freeStk, _AlgMatrixLLR::tol, and _AlgMatrixLLRE::val.

Referenced by AlgMatrixAdd(), AlgMatrixLLRSetAll(), AlgMatrixMul(), AlgMatrixScale(), AlgMatrixScaleAdd(), AlgMatrixSet(), AlgMatrixSub(), and AlgMatrixTranspose().

double AlgMatrixValue ( AlgMatrix  mat,
size_t  row,
size_t  col 
)

Returns the value in the matrix at the given coordinates.

Returns
Value in matrix at given coordinates.
Parameters
matGiven matrix.
rowGiven row.
colGiven column.

References ALG_MATRIX_LLR, ALG_MATRIX_RECT, ALG_MATRIX_SYM, AlgMatrixLLRValue(), _AlgMatrixSym::array, _AlgMatrix::core, _AlgMatrix::llr, _AlgMatrix::sym, and _AlgMatrixCore::type.

double AlgMatrixLLRValue ( AlgMatrixLLR mat,
size_t  row,
size_t  col 
)

Returns the value in the matrix at the given coordinates.

Returns
Value in matrix at given coordinates.
Parameters
matLinked list row matrix.
rowGiven row.
colGiven column.

References _AlgMatrixLLRE::col, _AlgMatrixLLRE::nxt, _AlgMatrixLLR::tbl, and _AlgMatrixLLRE::val.

Referenced by AlgMatrixMul(), AlgMatrixTrace(), and AlgMatrixValue().

AlgError AlgMatrixLLRExpand ( AlgMatrixLLR mat,
size_t  nE 
)

Ensures that there are at least the requested number of free linked list row matrix entries available.

Returns
Alg error code.
Parameters
matLinked list row matrix.
nERequested minimum number of free entries or if zero a default number of entries are added.

References AlcCalloc(), AlcFree(), AlcFreeStackPush(), ALG_ERR_MALLOC, ALG_ERR_NONE, AlgMatrixLLREFree(), _AlgMatrixLLR::blk, _AlgMatrixLLRE::col, _AlgMatrixLLR::freeStk, _AlgMatrixLLR::maxEnt, _AlgMatrixLLR::numEnt, _AlgMatrixLLRE::nxt, _AlgMatrixLLR::tbl, and _AlgMatrixLLRE::val.

Referenced by AlgMatrixLLRCopyInPlace(), AlgMatrixLLRNew(), AlgMatrixLLRSet(), AlgMatrixLLRSetAll(), AlgMatrixReadAscii(), and AlgMatrixScalar().

void AlgMatrixZero ( AlgMatrix  mat)
void AlgMatrixRectZero ( AlgMatrixRect mat)

Sets all elements of the matrix to zero.

Parameters
matGiven rectangular matrix.

References _AlgMatrixRect::array, _AlgMatrixRect::nC, and _AlgMatrixRect::nR.

Referenced by AlgMatrixScalar(), and AlgMatrixZero().

void AlgMatrixSymZero ( AlgMatrixSym mat)

Sets all elements of the matrix to zero.

Parameters
matGiven symmetric matrix.

References _AlgMatrixSym::array, and _AlgMatrixSym::nR.

Referenced by AlgMatrixScalar(), and AlgMatrixZero().

void AlgMatrixLLRZero ( AlgMatrixLLR mat)
AlgError AlgMatrixSetAll ( AlgMatrix  mat,
double  val 
)

Sets all elements of the given matrix to the given value.

Returns
Alg error code. An error may occur if the matrix type is not appropriate (ie not one of ALG_MATRIX_RECT, ALG_MATRIX_SYM or ALG_MATRIX_LLR) or if the elements of a ALG_MATRIX_LLR matrix can not be allocated.
Parameters
matGiven matrix.
valGiven value.

References ALG_ERR_MATRIX_TYPE, ALG_ERR_NONE, ALG_MATRIX_LLR, ALG_MATRIX_RECT, ALG_MATRIX_SYM, AlgMatrixLLRSetAll(), AlgMatrixRectSetAll(), AlgMatrixSymSetAll(), _AlgMatrix::core, _AlgMatrix::llr, _AlgMatrix::rect, _AlgMatrix::sym, and _AlgMatrixCore::type.

Referenced by WlzCMeshCompSurfMap().

void AlgMatrixRectSetAll ( AlgMatrixRect mat,
double  val 
)

Sets all elements of the given rectangular matrix to the given value.

Parameters
matGiven rectangular matrix.
valGiven value.

References _AlgMatrixRect::array, _AlgMatrixRect::nC, and _AlgMatrixRect::nR.

Referenced by AlgMatrixSetAll().

void AlgMatrixSymSetAll ( AlgMatrixSym mat,
double  val 
)

Sets all elements of the given symmetric matrix to the given value.

Parameters
matGiven symmetric matrix.
valGiven value.

References _AlgMatrixSym::array, and _AlgMatrixSym::nR.

Referenced by AlgMatrixSetAll().

AlgError AlgMatrixLLRSetAll ( AlgMatrixLLR mat,
double  val 
)

Sets all elements of the given linked list row matrix to the given value. This is not an efficient use of a linked list row matrix since all values will be allocate and set.

Returns
Error code.
Parameters
matGiven linked list row matrix.
valGiven value.

References ALG_ERR_NONE, AlgMatrixLLRExpand(), AlgMatrixLLRSet(), AlgMatrixLLRZero(), _AlgMatrixLLR::nC, _AlgMatrixLLR::nR, and _AlgMatrixLLR::tol.

Referenced by AlgMatrixSetAll().

void AlgMatrixLLRERemove ( AlgMatrixLLR mat,
size_t  row,
size_t  col 
)

Removes the entry at the given coordinates.

Parameters
matLinked list row matrix.
rowRow coordinate.
colColumn coordinate.

References AlgMatrixLLREFree(), _AlgMatrixLLRE::col, _AlgMatrixLLR::nR, _AlgMatrixLLRE::nxt, and _AlgMatrixLLR::tbl.

Referenced by AlgMatrixLLRSet(), and AlgMatrixScale().

AlgError AlgMatrixCGSolve ( AlgMatrix  aM,
double *  xV,
double *  bV,
AlgMatrix  wM,
void(*)(void *, AlgMatrix, double *, double *)  pFn,
void *  pDat,
double  tol,
int  maxItr,
double *  dstTol,
int *  dstItr 
)

Conjugate Gradient iterative method with preconditioning for the solution of linear systems with the form \(\mathbf{A} \mathbf{x} = \mathbf{b}\). \(\mathbf{A}\) must be a symmetric postive definite matrix, i.e. \({\mathbf{x}}^T \mathbf{A} \mathbf{x} < 0\), \(\forall \mathbf{x} \not= \mathbf{0}\), \(\mathbf{x} \in R^n\). Convergence is tested using: \( \frac{\| \mathbf{b} - mathbf{A} \mathbf{x} \|} {\|\mathbf{b}\|} < \delta\). If the preconditioning function pFn is non NULL then it is called, passing the preconditioning data pDat, as: (*pFn)(void *pDat, double **aM, double *r, double *z, int n) to solve \(\mathbf{A} \mathbf{z} = \mathbf{r}\) for \(\mathbf{z}\) with the solution overwriting the initial contents of z.

Returns
Error code.
Parameters
aMMatrix \(\mathbf{A}\).
xVMatrix \(\mathbf{x}\) which should contain an initial estimate although this may be \(\mathbf{0}\).
bVMatrix \(\mathbf{b}\).
wMMatrix with dimensions [4, n], this must be a rectangular matrix.
pFnPreconditioning function.
pDatData to be passed to preconditioning function.
tolTolerance required, \(\delta\).
maxItrThe maximum number of itterations.
dstTolDestination pointer for the residual after the final iteration, may be NULL.
dstItrDestination pointer for the actual number of itterations performed, may be NULL.

References ALG_ERR_CONVERGENCE, ALG_ERR_FUNC, ALG_ERR_NONE, ALG_MATRIX_LLR, ALG_MATRIX_RECT, ALG_MATRIX_SYM, AlgMatrixVectorMul(), AlgMatrixZero(), AlgVectorCopy(), AlgVectorDot(), AlgVectorNorm(), AlgVectorScaleAdd(), AlgVectorSub(), _AlgMatrixRect::array, _AlgMatrix::core, _AlgMatrixCore::nC, _AlgMatrixCore::nR, _AlgMatrix::rect, and _AlgMatrixCore::type.

Referenced by WlzAffineTransformLSqDQ3D(), WlzBasisFnScalarMOS3DFromCPts(), and WlzCMeshCompSurfMap().

AlgError AlgMatrixGaussSolve ( AlgMatrix  aMat,
double *  xMat 
)

Solves the matrix equation A.x = b for x by Gaussian elimination with partial pivoting. Matrix A is the matrix of coefficients.

Returns
Error code.
Parameters
aMatThe augmented matrix of size aSz x (aSz + 1), whose columns 0 - (aSz - 1) are the corresponding columns of A, and aSz'th column is b. Overwritten on exit.
xMatOn exit contains the solution matrix x.

References ALG_DBG, ALG_DBG_LVL_1, ALG_DBG_LVL_FN, ALG_ERR_FUNC, ALG_ERR_MATRIX_HOMOGENEOUS, ALG_ERR_MATRIX_SINGULAR, ALG_ERR_NONE, ALG_MATRIX_RECT, _AlgMatrixRect::array, _AlgMatrix::core, _AlgMatrixCore::nC, _AlgMatrixCore::nR, _AlgMatrixRect::nR, _AlgMatrix::rect, and _AlgMatrixCore::type.

Referenced by AlgPolynomialLSq().

AlgError AlgMatrixSolveLSQR ( AlgMatrix  aM,
double *  bV,
double *  xV,
double  damping,
double  relErrA,
double  relErrB,
long  maxItr,
long  condLim,
int *  dstTerm,
long *  dstItr,
double *  dstFNorm,
double *  dstCondN,
double *  dstResNorm,
double *  dstResNormA,
double *  dstNormX 
)

This function finds a solution x to the following problems:

Returns
    1. Unsymmetric equations --    solve  A*x = b

    2. Linear least squares  --    solve  A*x = b
                                   in the least-squares sense

    3. Damped least squares  --    solve  (   A    )*x = ( b )
                                          ( damp*I )     ( 0 )
                                   in the least-squares sense

    where 'A' is a matrix with 'm' rows and 'n' columns, 'b' is an
    'm'-vector, and 'damp' is a scalar.  (All quantities are real.)
    The matrix 'A' is intended to be large and sparse.  

    Notation
    --------
    The following quantities are used in discussing the subroutine
    parameters:

    'Abar'   =  (   A    ),          'bbar'  =  ( b )
                ( damp*I )                      ( 0 )

    'r'      =  b  -  A*x,           'rbar'  =  bbar  -  Abar*x

    'rnorm'  =  sqrt( norm(r)**2  +  damp**2 * norm(x)**2 )
             =  norm( rbar )

    'rel_prec'  =  the relative precision of floating-point arithmetic
                   on the machine being used.  Typically 2.22e-16
                   with 64-bit arithmetic.

    LSQR  minimizes the function 'rnorm' with respect to 'x'.

    References
    ----------
    C.C. Paige and M.A. Saunders,  LSQR: An algorithm for sparse
         linear equations and sparse least squares,
         ACM Transactions on Mathematical Software 8, 1 (March 1982),
         pp. 43-71.
    C.C. Paige and M.A. Saunders,  Algorithm 583, LSQR: Sparse
         linear equations and least-squares problems,
         ACM Transactions on Mathematical Software 8, 2 (June 1982),
         pp. 195-209.
    C.L. Lawson, R.J. Hanson, D.R. Kincaid and F.T. Krogh,
         Basic linear algebra subprograms for Fortran usage,
         ACM Transactions on Mathematical Software 5, 3 (Sept 1979),
         pp. 308-323 and 324-325.
Parameters
aMMatrix A with nR rows and nC columns. The values of A are not modified by this function.
bVVector b with nR entries which are modified by this function.
xVVector x for with initial guess at solution and for the return of the solution with nC entries.
dampingDamping parameter, set to 0.0 for no damping.
relErrAAn estimate of the relative error in matrix A. If A is accurate to about 6 digits, set to 1.0e-06.
relErrBAn estimate of the relative error in vector b. Set as for relErrA.
maxItrAn upper limit on the number of iterations. If zero set to nC.
condLimUpper limit on the condition number of matrix 'Abar'. Iterations will be terminated if a computed estimate of exceeds this condition number.
dstTermDestination pointer for termination code, which may be NULL. The values used are:
  • 0 x = x0 is the exact solution. No iterations were performed.
  • 1 The equations A*x = b are probably compatible. Norm(A*x - b) is sufficiently small, given the values of relErrA and relErrB.
  • 2 The system A*x = b is probably not compatible. A least-squares solution has been obtained that is sufficiently accurate, given the value of relErrA.
  • 3 An estimate of cond('Abar') has exceeded condLim. The system A*x = b appears to be ill-conditioned.
  • 4 The equations A*x = b are probably compatible. Norm(A*x - b) is as small as seems reasonable on this machine.
  • 5 The system A*x = b is probably not compatible. A least-squares solution has been obtained that is as accurate as seems reasonable on this machine.
  • 6 Cond('Abar') seems to be so large that there is no point in doing further iterations, given the precision of this machine.
  • 7 The iteration limit maxItr was reached.
dstItrDestination pointer for the number of iterations performed, which may be NULL.
dstFNormDestination pointer for an estimate of the Frobenius norm of 'Abar', may be NULL. This is the square-root of the sum of squares of the elements of 'Abar'. If 'damping' is small and if the columns of 'A' have all been scaled to have length 1.0, then the Frobenius norm should increase to roughly sqrt('nC').
dstCondNDestination pointer for an estimate of the condition number of 'Abar', may be NULL.
dstResNormDestination pointer for an estimate of the final value of norm('rbar'), the function being minimized. This will be small if A*x = b has a solution. May be NULL.
dstResNormADestination pointer for an estimate of the final value of the residual for the usual normal equations. This should be small in all cases. May be NULL.
dstNormXDestination pointer for an estimate of the final solution vector 'x'.

References AlcCalloc(), AlcFree(), ALG_ERR_CONVERGENCE, ALG_ERR_MALLOC, ALG_ERR_MATRIX_CONDITION, ALG_ERR_NONE, ALG_SQR, AlgMatrixTVectorMulAdd(), AlgMatrixVectorMulAdd(), AlgVectorCopy(), AlgVectorNorm(), AlgVectorScale(), _AlgMatrix::core, _AlgMatrixCore::nC, and _AlgMatrixCore::nR.

Referenced by WlzCMeshCompSurfMap().

AlgError AlgMatrixLUSolveRaw3 ( double **  aM,
double *  bV,
int  bSz 
)

Solves the matrix equation A.x = b for x, where A is a 3x3 libAlc double array matrix. On return the matrix A is overwritten with its LU decomposition and the column vector b is overwritten with the solution column matrix x.

Returns
Error code.
Parameters
aMRaw 3x3 array matrix A.
bVColumn vector b.
bSzSize (number of columns) of vector b (and x).

References ALG_ERR_NONE, AlgMatrixLUBackSubRaw(), and AlgMatrixLUDecompRaw().

Referenced by AlgMatrixLUInvertRaw3(), and Wlz3DViewGetIntersectionPoint().

AlgError AlgMatrixLUSolveRaw4 ( double **  aM,
double *  bV,
int  bSz 
)

Solves the matrix equation A.x = b for x, where A is a 4x4 libAlc double array matrix. On return the matrix A is overwritten with its LU decomposition and the column vector b is overwritten with the solution column matrix x.

Returns
Error code.
Parameters
aMRaw 4x4 array matrix A.
bVColumn vector b.
bSzSize (number of columns) of vector b (and x).

References ALG_ERR_NONE, AlgMatrixLUBackSubRaw(), and AlgMatrixLUDecompRaw().

Referenced by AlgMatrixLUInvertRaw4(), and WlzAffineTransformLSqRegWlz2D().

AlgError AlgMatrixLUSolve ( AlgMatrix  aM,
double *  bV,
int  bSz 
)

Solves the matrix equation A.x = b for x, where A is a square matrix. On return the matrix A is overwritten with its LU decomposition and the column vector b is overwritten with the solution column matrix x.

Returns
Error code.
Parameters
aMSquare matrix.
bVColumn vector b.
bSzSize (number of columns) of vector b (and x).

References ALG_ERR_FUNC, ALG_ERR_MATRIX_TYPE, ALG_ERR_NONE, ALG_MATRIX_RECT, AlgMatrixLUSolveRaw(), _AlgMatrixRect::array, _AlgMatrix::core, _AlgMatrixRect::nR, _AlgMatrix::rect, and _AlgMatrixCore::type.

AlgError AlgMatrixLUSolveRaw ( double **  aM,
int  aSz,
double *  bV,
int  bSz 
)

Solves the matrix equation A.x = b for x, where A is a square libAlc double array matrix. On return the matrix A is overwritten with its LU decomposition and the column vector b is overwritten with the solution column matrix x.

Returns
Error code.
Parameters
aMSquare libalc array matrix.
aSzSize of matrix A
bVColumn vector b.
bSzSize (number of columns) of vector b (and x).

References AlcFree(), AlcMalloc(), ALG_ERR_FUNC, ALG_ERR_MALLOC, ALG_ERR_NONE, AlgMatrixLUBackSubRaw(), and AlgMatrixLUDecompRaw().

Referenced by AlgMatrixLUInvertRaw(), and AlgMatrixLUSolve().

AlgError AlgMatrixLUInvertRaw3 ( double **  aM)

Calculates the inverse of a 3x3 libAlc double array matrix.

Returns
Error code.
Parameters
aMRaw double array of size 3x3.

References ALG_ERR_NONE, and AlgMatrixLUSolveRaw3().

AlgError AlgMatrixLUInvertRaw4 ( double **  aM)

Calculates the inverse of a 4x4 libAlc double array matrix.

Returns
Error code.
Parameters
aMRaw double array of size 4x4.

References ALG_ERR_NONE, and AlgMatrixLUSolveRaw4().

Referenced by WlzGeomTetraAffineSolveLU().

AlgError AlgMatrixLUInvert ( AlgMatrix  aM)

Calculates the inverse of a square matrix.

Returns
Error code.
Parameters
aMSquare libAlc double array matrix A.

References ALG_ERR_FUNC, ALG_ERR_MATRIX_TYPE, ALG_ERR_NONE, ALG_MATRIX_RECT, AlgMatrixLUInvertRaw(), _AlgMatrixRect::array, _AlgMatrix::core, _AlgMatrixRect::nR, _AlgMatrix::rect, and _AlgMatrixCore::type.

AlgError AlgMatrixLUInvertRaw ( double **  aM,
int  aSz 
)

Calculates the inverse of a square matrix.

Returns
Error code.
Parameters
aMSquare libAlc double array matrix A.
aSzSize of matrix A.

References AlcCalloc(), AlcFree(), ALG_ERR_FUNC, ALG_ERR_MALLOC, ALG_ERR_NONE, and AlgMatrixLUSolveRaw().

Referenced by AlgMatrixLUInvert(), and WlzAffineTransformInverse().

AlgError AlgMatrixLUDetermRaw3 ( double **  aM,
double *  det 
)

Calculates the determinant of a 3x3 double libAlc array matrix. The matrix is overwitten with its LU decomposition on return.

Returns
Error code.
Parameters
aMRaw 3x3 array matrix.
detDestination ptr for the determinant (may be NULL).

References ALG_ERR_NONE, and AlgMatrixLUDecompRaw().

AlgError AlgMatrixLUDetermRaw4 ( double **  aM,
double *  det 
)

Calculates the determinant of a 4x4 double libAlc array matrix. The matrix is overwitten with its LU decomposition on return.

Returns
Error code.
Parameters
aMRaw 4x4 array matrix.
detDestination ptr for the determinant (may be NULL).

References ALG_ERR_NONE, and AlgMatrixLUDecompRaw().

AlgError AlgMatrixLUDeterm ( AlgMatrix  aM,
double *  det 
)

Calculates the determinant of a square matrix. The matrix is overwitten with its LU decomposition on return.

Returns
Error code.
Parameters
aMMatrix A.
detDestination ptr for the determinant (may be NULL).

References ALG_ERR_FUNC, ALG_ERR_MATRIX_TYPE, ALG_ERR_NONE, ALG_MATRIX_RECT, AlgMatrixLUDetermRaw(), _AlgMatrixRect::array, _AlgMatrix::core, _AlgMatrixRect::nR, _AlgMatrix::rect, and _AlgMatrixCore::type.

AlgError AlgMatrixLUDetermRaw ( double **  aM,
int  aSz,
double *  det 
)

Calculates the determinant of a square double libAlc array matrix. The matrix is overwitten with its LU decomposition on return.

Returns
Error code.
Parameters
aMRaw matrix A.
aSzSize of matrix A.
detDestination ptr for the determinant (may be NULL).

References AlcFree(), AlcMalloc(), ALG_ERR_FUNC, ALG_ERR_MALLOC, ALG_ERR_NONE, and AlgMatrixLUDecompRaw().

Referenced by AlgMatrixLUDeterm(), and WlzCMeshStrainTensorAtPts().

AlgError AlgMatrixLUDecomp ( AlgMatrix  aM,
int *  iV,
double *  evenOdd 
)

Replaces the given matrix with the LU decomposition of a row-wise permutation of itself. The given index vector is used to record the permutation effected by the partial pivoting.

Returns
Error code.
Parameters
aMMatrix A.
iVIndex vector.
evenOddSet to +/-1.0 depending on whether the permutation was even or odd (may be NULL).

References ALG_ERR_FUNC, ALG_ERR_MATRIX_TYPE, ALG_ERR_NONE, ALG_MATRIX_RECT, AlgMatrixLUDecompRaw(), _AlgMatrixRect::array, _AlgMatrix::core, _AlgMatrixRect::nR, _AlgMatrix::rect, and _AlgMatrixCore::type.

AlgError AlgMatrixLUDecompRaw ( double **  aM,
int  aSz,
int *  iV,
double *  evenOdd 
)

Replaces the given matrix with the LU decomposition of a row-wise permutation of itself. The given index vector is used to record the permutation effected by the partial pivoting.

Returns
Error code.
Parameters
aMGiven raw matrix A.
aSzSize of matrix A.
iVIndex vector.
evenOddSet to +/-1.0 depending on whether the permutation was even or odd (may be NULL).

References AlcFree(), AlcMalloc(), ALG_ERR_FUNC, ALG_ERR_MALLOC, ALG_ERR_MATRIX_SINGULAR, and ALG_ERR_NONE.

Referenced by AlgMatrixLUDecomp(), AlgMatrixLUDetermRaw(), AlgMatrixLUDetermRaw3(), AlgMatrixLUDetermRaw4(), AlgMatrixLUSolveRaw(), AlgMatrixLUSolveRaw3(), AlgMatrixLUSolveRaw4(), WlzKrigOSetModelSV2D(), and WlzKrigOSetModelSV3D().

AlgError AlgMatrixLUBackSub ( AlgMatrix  aM,
int *  iV,
double *  bV 
)

Solves the set of of linear equations A.x = b where A is input as its LU decomposition determined with AlgMatrixLUDecomp() of a row-wise permutation of itself. The given index vector is used to record the permutation effected by the partial pivoting.

Returns
Error code.
Parameters
aMMatrix A.
iVIndex vector.
bVColumn vector b.

References ALG_ERR_FUNC, ALG_ERR_MATRIX_TYPE, ALG_ERR_NONE, ALG_MATRIX_RECT, AlgMatrixLUBackSubRaw(), _AlgMatrixRect::array, _AlgMatrix::core, _AlgMatrixRect::nR, _AlgMatrix::rect, and _AlgMatrixCore::type.

Referenced by WlzKrigOWeightsSolve().

AlgError AlgMatrixLUBackSubRaw ( double **  aM,
int  aSz,
int *  iV,
double *  bV 
)

Solves the set of of linear equations A.x = b where A is input as its LU decomposition determined with AlgMatrixLUDecompRaw() of a row-wise permutation of itself. The given index vector is used to record the permutation effected by the partial pivoting.

Returns
Error code.
Parameters
aMGiven raw libAlc array matrix A.
aSzSize of matrix A.
iVIndex vector.
bVColumn vector b.

References ALG_ERR_FUNC, and ALG_ERR_NONE.

Referenced by AlgMatrixLUBackSub(), AlgMatrixLUSolveRaw(), AlgMatrixLUSolveRaw3(), and AlgMatrixLUSolveRaw4().

void AlgMatrixAdd ( AlgMatrix  aM,
AlgMatrix  bM,
AlgMatrix  cM 
)

Computes the sum of two matricies and returns the result in a third supplied matrix:

\[ \mathbf{A} = \mathbf{B} + \mathbf{C} \]

The dimensions of the all three matricies must be nR, nC. It is safe to supply the same matrix as any combination of aM, bM and cM.

Returns
void
Note
For efficiency the given parameters are not checked.
Matrix size is limited only by address space.
All matrices must be of the same type (not checked for).
Parameters
aMSupplied matrix for result, \(\mathbf{A}\).
bMFirst matrix in the sum, \(\mathbf{B}\).
cMSecond matrix in the sum, \(\mathbf{C}\).

References ALG_MATRIX_LLR, ALG_MATRIX_RECT, ALG_MATRIX_SYM, AlgMatrixLLRSet(), AlgMatrixLLRZero(), _AlgMatrixRect::array, _AlgMatrixSym::array, _AlgMatrixLLRE::col, _AlgMatrix::core, _AlgMatrix::llr, _AlgMatrixCore::nC, _AlgMatrixCore::nR, _AlgMatrixLLRE::nxt, _AlgMatrix::rect, _AlgMatrix::sym, _AlgMatrixLLR::tbl, _AlgMatrixCore::type, and _AlgMatrixLLRE::val.

void AlgMatrixSub ( AlgMatrix  aM,
AlgMatrix  bM,
AlgMatrix  cM 
)

Subtracts on matrix from another and returns the result in a third supplied matrix:

\[ \mathbf{A} = \mathbf{B} - \mathbf{C} \]

The dimensions of the all three matricies must be nR, nC. It is safe to supply the same matrix as any combination of aM, bM and cM.

Returns
void
Note
For efficiency the given parameters are not checked.
Matrix size is limited only by address space.
Parameters
aMSupplied matrix for result, \(\mathbf{A}\).
bMFirst matrix in the subtraction, \(\mathbf{B}\).
cMSecond matrix in the subtraction, \(\mathbf{C}\).

References ALG_MATRIX_LLR, ALG_MATRIX_RECT, ALG_MATRIX_SYM, AlgMatrixLLRSet(), AlgMatrixLLRZero(), _AlgMatrixRect::array, _AlgMatrixSym::array, _AlgMatrixLLRE::col, _AlgMatrix::core, _AlgMatrix::llr, _AlgMatrixCore::nC, _AlgMatrixCore::nR, _AlgMatrixLLRE::nxt, _AlgMatrix::rect, _AlgMatrix::sym, _AlgMatrixLLR::tbl, _AlgMatrixCore::type, and _AlgMatrixLLRE::val.

Referenced by WlzAffineTransformLSqDQ2D(), and WlzAffineTransformLSqDQ3D().

void AlgMatrixMul ( AlgMatrix  aM,
AlgMatrix  bM,
AlgMatrix  cM 
)

Computes the product of two matricies and returns the result in a third supplied matrix:

\[ \mathbf{A} = \mathbf{B} \mathbf{C} \]

The dimensions of the result matrix (aM) must be cR, bC. All the matrices must be valid and of the same type except in the case of multiplying tow symmetric matrices when the result matrix must be rectangular (because in general the product of two symmetric matrices is not symmetric).

Returns
void
Note
For efficiency the given parameters are not fully checked.
Matrix size is limited only by address space.
Parameters
aMSupplied matrix for result, \(\mathbf{A}\)
bMFirst matrix in the product, \(\mathbf{B}\)
cMSecond matrix in the product, \(\mathbf{C}\)

References ALG_ERR_MATRIX_TYPE, ALG_ERR_NONE, ALG_MATRIX_LLR, ALG_MATRIX_RECT, ALG_MATRIX_SYM, AlgMatrixLLRSet(), AlgMatrixLLRValue(), _AlgMatrixRect::array, _AlgMatrixLLRE::col, _AlgMatrix::core, _AlgMatrix::llr, _AlgMatrixCore::nC, _AlgMatrixCore::nR, _AlgMatrixLLRE::nxt, _AlgMatrix::rect, _AlgMatrixLLR::tbl, _AlgMatrixCore::type, and _AlgMatrixLLRE::val.

Referenced by WlzAffineTransformLSqDQ2D(), and WlzAffineTransformLSqDQ3D().

double AlgMatrixTrace ( AlgMatrix  aM)

Computes the trace of the given matrix.

Returns
The trace of the matrix.
Note
For efficiency the given parameters are not checked.
Matrix size is limited only by address space.
Parameters
aMSupplied matrix.

References ALG_MATRIX_LLR, ALG_MATRIX_RECT, ALG_MATRIX_SYM, ALG_MIN, AlgMatrixLLRValue(), _AlgMatrixRect::array, _AlgMatrixSym::array, _AlgMatrix::core, _AlgMatrix::llr, _AlgMatrixCore::nC, _AlgMatrixCore::nR, _AlgMatrix::rect, _AlgMatrix::sym, and _AlgMatrixCore::type.

void AlgMatrixTranspose ( AlgMatrix  aM,
AlgMatrix  bM 
)

Computes the transpose of the given matrix:

\[ \mathbf{A} = \mathbf{B^T} \]

aM = transpose(bM). The dimensions of the result matrix (aM) must be aR == bC, aC == bR.

Returns
void
Note
For efficiency the given parameters are not checked.
Matrix size is limited only by address space.
Parameters
aMSupplied matrix for result, \(\mathbf{A}\).
bMMatrix to transpose, \(\mathbf{B}\).

References ALG_MATRIX_LLR, ALG_MATRIX_RECT, ALG_MATRIX_SYM, AlgMatrixLLRSet(), AlgMatrixLLRZero(), _AlgMatrixRect::array, _AlgMatrixSym::array, _AlgMatrixLLRE::col, _AlgMatrix::core, _AlgMatrix::llr, _AlgMatrixCore::nC, _AlgMatrixCore::nR, _AlgMatrixLLRE::nxt, _AlgMatrix::rect, _AlgMatrix::sym, _AlgMatrixLLR::tbl, _AlgMatrixCore::type, and _AlgMatrixLLRE::val.

Referenced by WlzAffineTransformLSqDQ2D(), and WlzAffineTransformLSqDQ3D().

void AlgMatrixCopy ( AlgMatrix  aM,
AlgMatrix  bM 
)

Copies the values of the matrix bM to the result matric aM:

\[ \mathbf{A} = \mathbf{B} \]

.

Returns
void
Note
For efficiency the given parameters are not checked.
Matrix size is limited only by address space.
Parameters
aMSupplied matrix for result, \(\mathbf{A}\).
bMMatrix to copy, \(\mathbf{B}\).

References ALG_MATRIX_LLR, ALG_MATRIX_RECT, ALG_MATRIX_SYM, AlgMatrixLLRCopyInPlace(), AlgVectorCopy(), _AlgMatrixRect::array, _AlgMatrixSym::array, _AlgMatrix::core, _AlgMatrix::llr, _AlgMatrixCore::nC, _AlgMatrixCore::nR, _AlgMatrix::rect, _AlgMatrix::sym, and _AlgMatrixCore::type.

Referenced by WlzAffineTransformLSqDQ3D().

void AlgMatrixScale ( AlgMatrix  aM,
AlgMatrix  bM,
double  sv 
)

Multiplies the given matrix by the given scalar:

\[ \mathbf{A} = s \mathbf{B} \]

.

Note
For efficiency the given parameters are not checked.
Matrix size is limited only by address space.
Parameters
aMSupplied matrix for result, \(\mathbf{A}\).
bMGiven matrix to scale, \(\mathbf{B}\).
svScalar value, \(s\).

References ALG_MATRIX_LLR, ALG_MATRIX_RECT, ALG_MATRIX_SYM, AlgMatrixLLRERemove(), AlgMatrixLLRSet(), AlgVectorScale(), _AlgMatrixRect::array, _AlgMatrixSym::array, _AlgMatrixLLRE::col, _AlgMatrix::core, _AlgMatrix::llr, _AlgMatrixCore::nC, _AlgMatrixCore::nR, _AlgMatrixLLR::nR, _AlgMatrixLLRE::nxt, _AlgMatrix::rect, _AlgMatrix::sym, _AlgMatrixLLR::tbl, _AlgMatrixLLR::tol, _AlgMatrixCore::type, and _AlgMatrixLLRE::val.

Referenced by WlzAffineTransformLSqDQ2D(), and WlzAffineTransformLSqDQ3D().

void AlgMatrixScaleAdd ( AlgMatrix  aM,
AlgMatrix  bM,
AlgMatrix  cM,
double  sv 
)

Multiplies the a matrix by a scalar and then adds another matrix:

\[ \mathbf{A} = \mathbf{B} + s \mathbf{C}. \]

All the matrices must have the same dimensions.

Returns
void
Note
For efficiency the given parameters are not checked.
Matrix size is limited only by address space.
Parameters
aMSupplied matrix for result, \(\mathbf{A}\).
bMGiven matrix to scale, \(\mathbf{B}\).
cMMatrix too add, \(\mathbf{C}\).
svScalar value, \(s\).

References ALG_MATRIX_LLR, ALG_MATRIX_RECT, ALG_MATRIX_SYM, AlgMatrixLLRSet(), AlgMatrixLLRZero(), _AlgMatrixRect::array, _AlgMatrixSym::array, _AlgMatrixLLRE::col, _AlgMatrix::core, _AlgMatrix::llr, _AlgMatrixCore::nC, _AlgMatrixCore::nR, _AlgMatrixLLRE::nxt, _AlgMatrix::rect, _AlgMatrix::sym, _AlgMatrixLLR::tbl, _AlgMatrixCore::type, and _AlgMatrixLLRE::val.

void AlgMatrixScalar ( AlgMatrix  aM,
double  sv 
)

Sets the elements of the given square matrix so that it is a scalar matrix:

\[ \mathbf{A} = s \mathbf{I} \]

.

Note
For efficiency the given parameters are not checked.
This function assumes that the matrix has been allocated by AlcDouble2Malloc().
Parameters
aMSupplied matrix for result, \(\mathbf{A}\).
svScalar value, \(s\).

References ALG_ERR_NONE, ALG_MATRIX_LLR, ALG_MATRIX_RECT, ALG_MATRIX_SYM, ALG_MIN, AlgMatrixLLRENew(), AlgMatrixLLRExpand(), AlgMatrixLLRZero(), AlgMatrixRectZero(), AlgMatrixSymZero(), _AlgMatrixRect::array, _AlgMatrixLLRE::col, _AlgMatrix::core, _AlgMatrix::llr, _AlgMatrixCore::nC, _AlgMatrixCore::nR, _AlgMatrixLLRE::nxt, _AlgMatrix::rect, _AlgMatrix::sym, _AlgMatrixLLR::tbl, _AlgMatrixCore::type, and _AlgMatrixLLRE::val.

void AlgMatrixVectorMul ( double *  aV,
AlgMatrix  bM,
double *  cV 
)

Multiplies the matrix \(\mathbf{B}\) by the vector \(\mathbf{c}\):

\[ \mathbf{a} = \mathbf{B} \mathbf{c} \]

.

Note
This function assumes that the matrix has been allocated by either AlcDouble2Malloc() or AlcSymDouble2Malloc().
Matrix size is limited only by address space.
Parameters
aVSupplied vector for result.
bMMatrix \(mathbf{B}\).
cVVector \(\mathbf{c}\). \(mathbf{B}\).

References ALG_MATRIX_LLR, ALG_MATRIX_RECT, ALG_MATRIX_SYM, _AlgMatrixRect::array, _AlgMatrixSym::array, _AlgMatrixLLRE::col, _AlgMatrix::core, _AlgMatrix::llr, _AlgMatrixRect::nC, _AlgMatrixRect::nR, _AlgMatrixSym::nR, _AlgMatrixLLR::nR, _AlgMatrixLLRE::nxt, _AlgMatrix::rect, _AlgMatrix::sym, _AlgMatrixLLR::tbl, _AlgMatrixCore::type, and _AlgMatrixLLRE::val.

Referenced by AlgMatrixCGSolve(), and WlzCMeshCompSurfMap().

void AlgMatrixVectorMulAdd ( double *  aV,
AlgMatrix  bM,
double *  cV,
double *  dV 
)

Multiplies the matrix \(\mathbf{B}\) by the vector \(\mathbf{c}\) and adds the vector \(\mathbf{d}\):

\[ \mathbf{a} = \mathbf{B} \mathbf{c} + \mathbf{d} \]

.

Note
This function assumes that the matrix has been allocated by either AlcDouble2Malloc() or AlcSymDouble2Malloc().
Matrix size is limited only by address space.
Parameters
aVSupplied vector for result.
bMMatrix \(mathbf{B}\).
cVVector \(\mathbf{c}\).
dVVector \(\mathbf{d}\). \(mathbf{B}\).

References ALG_MATRIX_LLR, ALG_MATRIX_RECT, ALG_MATRIX_SYM, _AlgMatrixRect::array, _AlgMatrixSym::array, _AlgMatrixLLRE::col, _AlgMatrix::core, _AlgMatrix::llr, _AlgMatrixRect::nC, _AlgMatrixRect::nR, _AlgMatrixSym::nR, _AlgMatrixLLR::nR, _AlgMatrixLLRE::nxt, _AlgMatrix::rect, _AlgMatrix::sym, _AlgMatrixLLR::tbl, _AlgMatrixCore::type, and _AlgMatrixLLRE::val.

Referenced by AlgMatrixSolveLSQR().

void AlgMatrixVectorMulWAdd ( double *  aV,
AlgMatrix  bM,
double *  cV,
double *  dV,
double  s,
double  t 
)

Multiplies the matrix \(\mathbf{B}\) by the vector \(\mathbf{c}\) and adds the vector \(\mathbf{d}\) using the given weights \(s\) and \(t\):

\[ \mathbf{a} = s \mathbf{B} \mathbf{c} + t \mathbf{d} \]

.

Note
This function assumes that the matrix has been allocated by either AlcDouble2Malloc() or AlcSymDouble2Malloc().
Matrix size is limited only by address space.
Parameters
aVSupplied vector for result.
bMMatrix \(mathbf{B}\).
cVVector \(\mathbf{c}\).
dVVector \(\mathbf{d}\).
sFirst weighting scalar \(s\).
tSecond weighting scalar \(t\).

References ALG_MATRIX_LLR, ALG_MATRIX_RECT, ALG_MATRIX_SYM, _AlgMatrixRect::array, _AlgMatrixSym::array, _AlgMatrixLLRE::col, _AlgMatrix::core, _AlgMatrix::llr, _AlgMatrixRect::nC, _AlgMatrixRect::nR, _AlgMatrixSym::nR, _AlgMatrixLLR::nR, _AlgMatrixLLRE::nxt, _AlgMatrix::rect, _AlgMatrix::sym, _AlgMatrixLLR::tbl, _AlgMatrixCore::type, and _AlgMatrixLLRE::val.

void AlgMatrixTVectorMul ( double *  aV,
AlgMatrix  bM,
double *  cV 
)

Multiplies the transpose of matrix \(\mathbf{B}\) by the vector \(\mathbf{c}\):

\[ \mathbf{a} = \mathbf{B}^T \mathbf{c} \]

.

Note
This function assumes that the matrix has been allocated by either AlcDouble2Malloc() or AlcSymDouble2Malloc().
Matrix size is limited only by address space.
Parameters
aVSupplied vector for result.
bMMatrix \(mathbf{B}\).
cVVector \(\mathbf{c}\). \(mathbf{B}\).

References ALG_MATRIX_LLR, ALG_MATRIX_RECT, ALG_MATRIX_SYM, AlgVectorZero(), _AlgMatrixRect::array, _AlgMatrixSym::array, _AlgMatrixLLRE::col, _AlgMatrix::core, _AlgMatrix::llr, _AlgMatrixRect::nC, _AlgMatrixLLR::nC, _AlgMatrixRect::nR, _AlgMatrixSym::nR, _AlgMatrixLLR::nR, _AlgMatrixLLRE::nxt, _AlgMatrix::rect, _AlgMatrix::sym, _AlgMatrixLLR::tbl, _AlgMatrixCore::type, and _AlgMatrixLLRE::val.

void AlgMatrixTVectorMulAdd ( double *  aV,
AlgMatrix  bM,
double *  cV,
double *  dV 
)

Multiplies the transpose of matrix \(\mathbf{B}\) by the vector \(\mathbf{c}\):

\[ \mathbf{a} = \mathbf{B}^T \mathbf{c} + \mathbf{d} \]

.

Note
This function assumes that the matrix has been allocated by either AlcDouble2Malloc() or AlcSymDouble2Malloc().
Matrix size is limited only by address space.
Parameters
aVSupplied vector for result.
bMMatrix \(mathbf{B}\).
cVVector \(\mathbf{c}\).
dVVector \(\mathbf{d}\). \(mathbf{B}\).

References ALG_MATRIX_LLR, ALG_MATRIX_RECT, ALG_MATRIX_SYM, AlgVectorCopy(), _AlgMatrixSym::array, _AlgMatrixLLRE::col, _AlgMatrix::core, _AlgMatrix::llr, _AlgMatrixRect::nC, _AlgMatrixLLR::nC, _AlgMatrixRect::nR, _AlgMatrixSym::nR, _AlgMatrixLLR::nR, _AlgMatrixLLRE::nxt, _AlgMatrix::rect, _AlgMatrix::sym, _AlgMatrixLLR::tbl, _AlgMatrixCore::type, and _AlgMatrixLLRE::val.

Referenced by AlgMatrixSolveLSQR().

AlgError AlgMatrixRawInv2x2 ( double *  a00,
double *  a01,
double *  a10,
double *  a11 
)

Inverts a raw 2 x 2 matrix. All matrix values are passed using pointers and their values are set to those of the inverse matrix on return. If the given matrix is singular the ALG_ERR_MATRIX_SINGULAR error code is returned and no values are changed.

Returns
Alg error code, either ALG_ERR_NONE or ALG_ERR_MATRIX_SINGULAR.
Parameters
a00Matrix element row 0, column 0.
a01Matrix element row 0, column 1.
a10Matrix element row 1, column 0.
a11Matrix element row 1, column 1.

References ALG_ERR_MATRIX_SINGULAR, and ALG_ERR_NONE.

AlgError AlgMatrixRawInv3x3 ( double *  a00,
double *  a01,
double *  a02,
double *  a10,
double *  a11,
double *  a12,
double *  a20,
double *  a21,
double *  a22 
)

Inverts a raw 2 x 2 matrix. All matrix values are passed using pointers and their values are set to those of the inverse matrix on return. If the given matrix is singular the ALG_ERR_MATRIX_SINGULAR error code is returned and no values are changed.

Returns
Alg error code, either ALG_ERR_NONE or ALG_ERR_MATRIX_SINGULAR.
Parameters
a00Matrix element row 0, column 0.
a01Matrix element row 0, column 1.
a02Matrix element row 0, column 2.
a10Matrix element row 1, column 0.
a11Matrix element row 1, column 1.
a12Matrix element row 1, column 2.
a20Matrix element row 2, column 0.
a21Matrix element row 2, column 1.
a22Matrix element row 2, column 2.

References ALG_ERR_MATRIX_SINGULAR, and ALG_ERR_NONE.

AlgError AlgMatrixRSEigen ( AlgMatrix  aM,
double *  vM,
int  reqEV 
)

Determines the eigenvalues and eigenvectors of a real symmetric matrix by calling AlgMatrixRSTDiag() to create a tridiagonal symmetric matrix and then AlgMatrixTDiagQLI() to compute its eigenvalues and eigenvectors. The eigenvectors and eigenvalues are returned in descending eigenvalue order. For efficiency, the eigenvectors should only be computed if required.

Returns
Error code.
Parameters
aMGiven real symmetric matrix which contains the eigenvectors in it's columns on return if required.
vMGiven vector for the return of the eigenvalues.
reqEVNon zero if the eigenvectors are required.

References AlcFree(), AlcMalloc(), ALG_ERR_FUNC, ALG_ERR_MALLOC, ALG_ERR_NONE, ALG_MATRIX_RECT, AlgMatrixRSTDiag(), AlgMatrixTDiagQLI(), _AlgMatrixRect::array, _AlgMatrix::core, _AlgMatrixCore::nC, _AlgMatrixCore::nR, _AlgMatrix::rect, and _AlgMatrixCore::type.

Referenced by WlzAffineTransformLSqDQ2D(), WlzAffineTransformLSqDQ3D(), WlzCMeshExpansion(), and WlzCMeshStrainTensorAtPts().

AlgError AlgMatrixRSTDiag ( AlgMatrix  aMat,
double *  dVec,
double *  oVec 
)

An implementation of Householder's alorithm which reduces a real aSz x aSz symmetric matrix to symmetric tridiagonal form by aSz - 2 orthogonal similarity transformations. This is a translation of the functions/subroutines 'tred2' in the fortran EISPACK library and the Numerical Recipies book, both of these are descended from the algol procedure tred2 in Wilkinson J.H and Reinsch C. Linear Algebra, Volume II Handbook for Automatic Computation. Springer-Verlag, 1971. See Numerical Recipies in C: The Art of Scientific Computing. Cambridge University Press, 1992.

Returns
Error code.
Parameters
aMatGiven real symmetric matrix.
dVecGiven array for the return of the diagonal elements.
oVecGiven array for the return of the off diagonal elements with the first element set to 0.0.

References ALG_ERR_FUNC, ALG_ERR_NONE, ALG_MATRIX_RECT, _AlgMatrixRect::array, _AlgMatrix::core, _AlgMatrixCore::nC, _AlgMatrixCore::nR, _AlgMatrixRect::nR, _AlgMatrix::rect, and _AlgMatrixCore::type.

Referenced by AlgMatrixRSEigen().

AlgError AlgMatrixSVSolve ( AlgMatrix  aMat,
double *  bVec,
double  tol,
int *  dstIC 
)

Solves the matrix equation A.x = b for x, where A is a matrix with at least as many rows as columns. On return the matrix A is overwritten by the matrix U in the singular value decomposition: A = U.W.V'.

Returns
Error code.
Parameters
aMatMatrix A.
bVecColumn matrix b, overwritten by matrix x on return.
tolTolerance for singular values, 1.0e-06 should be suitable as a default value.
dstICDestination pointer for ill-conditioned flag, which may be NULL, but if not NULL is set to the number of singular values which are smaller than the maximum singular value times the given threshold.

References AlcCalloc(), AlcFree(), ALG_DBG, ALG_DBG_LVL_1, ALG_DBG_LVL_FN, ALG_ERR_FUNC, ALG_ERR_MALLOC, ALG_ERR_NONE, ALG_MATRIX_RECT, AlgMatrixRectFree(), AlgMatrixRectNew(), AlgMatrixSVBackSub(), AlgMatrixSVDecomp(), _AlgMatrix::core, _AlgMatrixCore::nC, _AlgMatrixRect::nC, _AlgMatrixCore::nR, _AlgMatrix::rect, and _AlgMatrixCore::type.

Referenced by WlzAffineTransformLSqDQ3D(), and WlzGeomCurvature().

AlgError AlgMatrixSVDecomp ( AlgMatrix  aMat,
double *  wVec,
AlgMatrix  vMat 
)

Performs singular value decomposition of the given matrix ( \(A\) ) and computes two additional matricies ( \(U\) ) and ( \(V\) ) such that:

\[ A = U W V^T \]

where \(V^T\) is the transpose of \(V\). The code for AlgMatrixSVDecomp was derived from: Numerical Recipies function svdcmp(), EISPACK subroutine SVD(), CERN subroutine SVD() and ACM algorithm 358. The singular values ( \(W\) ) and singular vectors ( \(V\) ) are returned unsorted. See AlgMatrixSVSolve() for a usage example.

Returns
Error code.
Parameters
aMatThe given matrix A, and U on return.
wVecThe diagonal matrix of singular values, returned as a vector.
vMatThe matrix V (not it's transpose).

References AlcCalloc(), AlcFree(), ALG_DBG, ALG_DBG_LVL_1, ALG_DBG_LVL_FN, ALG_ERR_FUNC, ALG_ERR_MALLOC, ALG_ERR_MATRIX_SINGULAR, ALG_ERR_NONE, ALG_MATRIX_RECT, _AlgMatrixRect::array, _AlgMatrix::core, _AlgMatrixCore::nC, _AlgMatrixRect::nC, _AlgMatrixCore::nR, _AlgMatrixRect::nR, _AlgMatrix::rect, and _AlgMatrixCore::type.

Referenced by AlgMatrixSVSolve(), WlzAffineTransformLSqReg2D(), WlzAffineTransformLSqReg3D(), WlzBasisFnConf2DFromCPts(), WlzBasisFnGauss2DFromCPts(), WlzBasisFnIMQ2DFromCPts(), WlzBasisFnIMQ3DFromCPts(), WlzBasisFnMQ2DFromCPts(), WlzBasisFnMQ3DFromCPts(), WlzBasisFnPoly2DFromCPts(), WlzBasisFnScalarMOS3DFromCPts(), WlzBasisFnTPS2DFromCPts(), WlzCMeshExpValues(), WlzFitPlaneSVD(), and WlzGeometryLSqOPlane().

AlgError AlgMatrixSVBackSub ( AlgMatrix  uMat,
double *  wVec,
AlgMatrix  vMat,
double *  bVec 
)

Solves the set of of linear equations A.x = b where A is input as its singular value decomposition in the three matricies U, W and V, as returned by AlgMatrixSVDecomp(). The code for AlgMatrixSVBackSub was derived from: Numerical Recipies function svbksb().

Returns
Error code.
Parameters
uMatGiven matrix U with nM rows and nN columns..
wVecThe diagonal matrix of singular values, returned as a vector with nN elements.
vMatThe matrix V (not it's transpose) with nN columns.
bVecColumn matrix b with nM elements, overwritten by column matrix x on return.

References AlcCalloc(), AlcFree(), ALG_DBG, ALG_DBG_LVL_1, ALG_DBG_LVL_FN, ALG_ERR_FUNC, ALG_ERR_MALLOC, ALG_ERR_NONE, ALG_MATRIX_RECT, _AlgMatrixRect::array, _AlgMatrix::core, _AlgMatrixCore::nC, _AlgMatrixRect::nC, _AlgMatrixCore::nR, _AlgMatrixRect::nR, _AlgMatrix::rect, and _AlgMatrixCore::type.

Referenced by AlgMatrixSVSolve(), WlzBasisFnConf2DFromCPts(), WlzBasisFnGauss2DFromCPts(), WlzBasisFnIMQ2DFromCPts(), WlzBasisFnIMQ3DFromCPts(), WlzBasisFnMQ2DFromCPts(), WlzBasisFnMQ3DFromCPts(), WlzBasisFnPoly2DFromCPts(), WlzBasisFnScalarMOS3DFromCPts(), and WlzBasisFnTPS2DFromCPts().

AlgError AlgMatrixTDiagQLI ( double *  dM,
double *  oM,
int  aSz,
AlgMatrix  zM 
)

Determines the eigenvalues and eigenvectors of a real symmetric tridiagonal matrix using the QL algorithm with implicit shifts. This is a based on the function 'tqli' in Numerical Recipies in C: The Art of Scientific Computiung. Cambridge University Press, 1992.

Returns
Error code.
Parameters
dMGiven diagonal elements of the tridiagonal matrix. On return it contains the eigenvalues.
oMGiven off diagonal elements of the tridiagonal matrix, with an arbitrary first element. On return it's contents are destroyed.
aSzSize of the array.
zMGiven array for the return of the eigenvectors, may be NULL if the eigenvectors are not required. If required the i'th eigenvector is returned in the i'th column of zM.

References ALG_ERR_CONVERGENCE, ALG_ERR_FUNC, ALG_ERR_NONE, ALG_MATRIX_RECT, _AlgMatrixRect::array, _AlgMatrix::core, _AlgMatrix::rect, and _AlgMatrixCore::type.

Referenced by AlgMatrixRSEigen().