BzrCrvInterp

(bzrintrp.c:211)

Prototype:

  void BzrCrvInterp(IrtRType *Result, IrtRType *Input, int Size)


Description:

Blends the input points using the Interp array and puts the result blended point in Result array. Input and Result array are of size Size.

Parameters:

Result: Where the interpolated control points will be placed.
Input: Points to interpolate at node parameter values.
Size: Of control polygon.


Returned Value:

void


Keywords:

interpolation


JacobiMatrixDiag4x4

(diag_mat.c:110)

Prototype:

  void JacobiMatrixDiag4x4(IrtRType M[4][4],
                           IrtRType U[4][4],
                           IrtRType D[4][4],
                           IrtRType V[4][4])


Description:

performs a diagonalization of a 4x4 matrix, as U D V, D diagonal and U and V orthogonal matrices U = V^-1 = V^T.

Parameters:

M: Matrix to diagonalize.
U: The left orthogonal (=rotation) matrix.
D: The diagonal matrix.
V: The eight orthogonal matrix.


Returned Value:

void


See Also:

JacobiMatrixDiagNxN

Keywords:




JacobiMatrixDiagNxN

(diag_mat.c:51)

Prototype:

  void JacobiMatrixDiagNxN(IrtRType *M[],
                           IrtRType *U[],
                           IrtRType *D[],
                           IrtRType *V[],
                           int n)


Description:

Performs a diagonalization of a NxN matrix, as U D V, D diagonal and U and V orthogonal matrices U = V^-1 = V^T.

Parameters:

M: Matrix to diagonalize.
U: The left orthogonal (=rotation) matrix.
D: The diagonal matrix.
V: The eight orthogonal matrix.
n: Dimension of the matrices as n times n.


Returned Value:

void


See Also:

JacobiMatrixDiag4x4

Keywords:




SvdDecomp

(nure_svd.c:342)

Prototype:

  void SvdDecomp(IrtRType **a, int m, int n, IrtRType *w, IrtRType **v)


Description:

SVD decomposition of matrix A as U W V, U of size m x m, V of size n x n and W is diagonal with the signular values, stored as a vector. Due to the Fortran style of this code all indices are from 1 to k, so an array from 1 to k actually have k + 1 elements in this C translation.

Parameters:

a: Input matrix also serving to store U orthogonal output, of dimension m x n, m >= n.
m, n: The dimensions of a.
w: The diagonal singular values, as a vector of length m.
v: The V orthogonal output, dimensions n x n.


Returned Value:

void


See Also:

SvdMatrixNxN

Keywords:




SvdLeastSqr

(nure_svd.c:798)

Prototype:

  IrtRType SvdLeastSqr(XtraSvdLSGenInfoStruct *SLSGI, IrtRType *x, IrtRType *b)


Description:

Least square solves A x = b. The vector X is of size Nx, vector b is of size NData and matrix A (SVD decomposed in SLSGI) is of size Nx by NData. Uses singular value decomposition.

Parameters:

SLSGI: The SVD least sqaure general info reference.
x: The vector of sought solution of size Nx.
b: The vector of coefficients of size NData.


Returned Value:

IrtRType: The reciprocal of the condition number, if A != NULL, zero otherwise.


See Also:

SvdMatrix4x4 SvdLeastSqrInit SvdLeastSqrFree

Keywords:

singular value decomposition linear systems


SvdLeastSqrFree

(nure_svd.c:845)

Prototype:

  void SvdLeastSqrFree(XtraSvdLSGenInfoStruct *SLSGI)


Description:

Least square solves A x = b. The vector X is of size Nx, vector b is of size NData and matrix A (SVD decomposed in SLSGI) is of size Nx by NData. Uses singular value decomposition.

Parameters:

SLSGI: The SVD least sqaure general info reference.


Returned Value:

void


See Also:

SvdMatrix4x4 SvdLeastSqr SvdLeastSqrInit

Keywords:

singular value decomposition linear systems


SvdLeastSqrInit

(nure_svd.c:658)

Prototype:

  XtraSvdLSGenInfoStruct *SvdLeastSqrInit(IrtRType *A,
                                          int NData,
                                          int Nx,
                                          IrtRType *RecipCondNum)


Description:

Init the least squares solver of A x = b. Matrix A is of size Nx by NData. Uses singular value decomposition.

Parameters:

A: The matrix of size Nx by NData.
NData, Nx: Dimensions of input.
RecipCondNum: The Recirocating value of the condition number.


Returned Value:

XtraSvdLSGenInfoStruct *: The SVD least sqaure general info reference.


See Also:

SvdMatrix4x4. SvdLeastSqr SvdLeastSqrFree

Keywords:

singular value decomposition linear systems


SvdMatrix4x4

(nure_svd.c:45)

Prototype:

  int SvdMatrix4x4(IrtHmgnMatType M,
                   IrtRType U[3][3],
                   IrtVecType S,
                   IrtRType V[3][3])


Description:

Perform singular value decomposition on the 3x3 main minor of the affine transformation matrix.

Parameters:

M: Source matrix to decompose.
U: Left orthonormal matrix.
S: Singular components.
V: Right orthonormal matrix.


Returned Value:

int: FALSE if failed, TRUE otherwise.


See Also:

SvdLeastSqr SvdMatrixNxN

Keywords:

singular value decomposition linear systems matrices


SvdMatrixNxN

(nure_svd.c:284)

Prototype:

  void SvdMatrixNxN(IrtRType *M, IrtRType *U, IrtRType *S, IrtRType *V, int n)


Description:

Perform singular value decomposition on an NxN matrix M as U S V.

Parameters:

M: Source matrix to decompose.
U: Left orthonormal matrix, size n x n.
S: A vector holding the singular components.
V: Right orthonormal matrix, size n.
n: size of square matrix, size n x n.


Returned Value:

void


See Also:

SvdLeastSqr SvdDecomp SvdMatrix3x3

Keywords:

singular value decomposition linear systems matrices