Ravelin
|
Linear algebra routines. More...
#include <LinAlg.h>
Public Types | |
enum | SVD { eSVD1, eSVD2 } |
Public Member Functions | |
void | compress () |
Compresses all memory. | |
void | free_memory () |
Frees all allocated memory. | |
MATRIXN & | pseudo_invert (MATRIXN &A, REAL tol=(REAL)-1.0) |
void | update_QR_rank1 (MATRIXN &Q, MATRIXN &R, const VECTORN &u, const VECTORN &v) |
Updates a QR factorization by a rank-1 update. More... | |
void | update_QR_delete_cols (MATRIXN &Q, MATRIXN &R, unsigned k, unsigned p) |
Updates a QR factorization by deleting p columns starting at column idx k. More... | |
void | update_QR_insert_cols (MATRIXN &Q, MATRIXN &R, MATRIXN &U, unsigned k) |
Updates a QR factorization by inserting one or more columns at column idx k. More... | |
void | update_QR_insert_rows (MATRIXN &Q, MATRIXN &R, MATRIXN &U, unsigned k) |
Updates a QR factorization by inserting a block of rows, starting at index k. More... | |
void | update_QR_delete_rows (MATRIXN &Q, MATRIXN &R, unsigned k, unsigned p) |
Updates a QR factorization by deleting a block of rows. More... | |
template<class X > | |
X & | inverse_LU (X &M, const std::vector< int > &pivwork) |
Computes a matrix inverse using the factorization determined via factor_LU() More... | |
template<class X > | |
unsigned | calc_rank (X &A, REAL tol=(REAL)-1.0) |
Calculates the rank of a matrix. | |
template<class Y > | |
MATRIXN & | nullspace (Y &A, MATRIXN &nullspace, REAL tol=-1.0) |
Computes the nullspace of a matrix. More... | |
template<class X > | |
REAL | cond (X &A) |
Computes the condition number of a matrix. | |
template<class X > | |
bool | is_SPSD (X &m, REAL tol) |
Determines whether a symmetric matrix is positive semi-definite. | |
template<class X > | |
bool | is_SPD (X &m, REAL tol) |
Determines whether a matrix is positive-definite. | |
template<class X , class Y > | |
void | eig_symm (X &A, Y &evals) |
Computes the eigenvalues of the matrix A. More... | |
template<class X , class Y > | |
void | eig_symm_plus (X &A_evecs, Y &evals) |
Computes the eigenvalues and eigenvectors of the matrix A. More... | |
template<class X , class MatU , class VecS , class MatV > | |
void | svd (X &A, MatU &U, VecS &S, MatV &V) |
template<class X , class MatU , class VecS , class MatV > | |
void | svd1 (X &A, MatU &U, VecS &S, MatV &V) |
Does an 'in place' SVD (destroying A) More... | |
template<class X , class MatU , class VecS , class MatV > | |
void | svd2 (X &A, MatU &U, VecS &S, MatV &V) |
Does an 'in place' SVD (destroying A), using divide and conquer algorithm. More... | |
template<class X > | |
X & | solve_symmetric_fast (MATRIXN &A, X &XB) |
Solves a symmetric, indefinite square matrix. More... | |
template<class X > | |
X & | inverse_symmetric (X &A) |
Inverts the symmetric, indefinite matrix A. More... | |
template<class X > | |
X & | invert (X &A) |
Inverts the matrix A using LU factorization. More... | |
template<class X , class Y , class Vec , class Z > | |
X & | solve_LS_fast (const Y &U, const Vec &S, const Z &V, X &XB, REAL tol=(REAL)-1.0) |
Most robust system of linear equations solver (solves AX = B) More... | |
template<class X , class Y > | |
X & | solve_LS_fast (Y &A, X &XB, SVD svd_algo, REAL tol) |
Most robust system of linear equations solver (solves AX = B) More... | |
template<class X , class Y > | |
Y & | solve_fast (X &A, Y &XB) |
Solves the general system AX = B. More... | |
template<class X , class Y > | |
X & | solve_LS_fast1 (Y &A, X &XB, REAL tol=(REAL)-1.0) |
template<class X , class Y > | |
X & | solve_LS_fast2 (Y &A, X &XB, REAL tol=(REAL)-1.0) |
template<class ARMat , class QMat > | |
void | factor_QR (ARMat &AR, QMat &Q, std::vector< int > &PI) |
Performs the QR factorization of a matrix with column pivoting. More... | |
template<class ARMat , class QMat > | |
void | factor_QR (ARMat &AR, QMat &Q) |
Performs the QR factorization of a matrix. More... | |
Static Public Member Functions | |
static void | factor_LDL (MATRIXN &M, std::vector< int > &IPIV) |
static void | givens (REAL a, REAL b, REAL &c, REAL &s) |
Less robust least squares solver (solves Ax = b) More... | |
static MATRIX2 | givens (REAL c, REAL s) |
Computes the givens matrix given a c and s. | |
static void | householder (REAL alpha, const VECTORN &x, REAL &tau, VECTORN &v) |
Computes a householder vector. | |
static VECTORN & | solve_sparse_direct (const SPARSEMATRIXN &A, const VECTORN &b, Transposition trans, VECTORN &x) |
static MATRIXN & | solve_sparse_direct (const SPARSEMATRIXN &A, const MATRIXN &B, Transposition trans, MATRIXN &X) |
template<class X > | |
static X & | solve_tridiagonal_fast (VECTORN &dl, VECTORN &d, VECTORN &du, X &XB) |
Solves a tridiagonal system. More... | |
template<class X > | |
static bool | factor_chol (X &A) |
Performs the Cholesky factorization of a matrix. More... | |
template<class X > | |
static bool | factor_LU (X &A, std::vector< int > &pivwork) |
Performs the LU factorization of a matrix. More... | |
template<class X , class Y > | |
static X & | solve_tri_fast (Y &A, bool utri, bool transpose_A, X &XB) |
Solves a triangular system of linear equations. More... | |
template<class X > | |
static X & | solve_LDL_fast (const MATRIXN &M, const std::vector< int > &pivwork, X &XB) |
Solves systems of linear equations using the factorization determined via factor_LDL() More... | |
template<class X , class Y > | |
static X & | solve_chol_fast (const Y &M, X &XB) |
Solves a system of linear equations using the factorization determined via factor_chol() More... | |
template<class Y , class X > | |
static X & | solve_LU_fast (const Y &M, bool transpose, const std::vector< int > &pivwork, X &XB) |
Solves a system of linear equations using the factorization determined via factor_LU() More... | |
template<class Y , class X > | |
static X & | solve_SPD_fast (Y &A, X &XB) |
Solves a system of equations A*X = B using a symmetric, positive-definite square matrix. More... | |
template<class X > | |
static X & | inverse_chol (X &A) |
Inverts the symmetric, positive-definite matrix A using Cholesky factorization. More... | |
template<class X > | |
static X & | inverse_SPD (X &A) |
Inverts the symmetric, positive-definite matrix A using Cholesky factorization. More... | |
Public Attributes | |
FastThreadable< MATRIXN > | workM |
work matrices | |
FastThreadable< MATRIXN > | workM2 |
FastThreadable< MATRIXN > | U |
work matrix (for SVD) | |
FastThreadable< std::vector < INTEGER > > | pivwork |
work STL integer vector (pivoting) | |
FastThreadable< MATRIXN > | V |
work matrix (for SVD) | |
FastThreadable< VECTORN > | S |
work vector (for SVD/eigenvalues) | |
FastThreadable< VECTORN > | workv |
work vectors (LAPACK routines) | |
FastThreadable< VECTORN > | workv2 |
FastThreadable< std::vector < INTEGER > > | iworkv |
work STL integer vector (LAPACK routines) | |
Linear algebra routines.
LinAlg is a set of static routines that interface to LAPACK. I have included only very few routines here, however they should be some of the most utilized: SVD, (SVD-based) pseudo-inverse, linear equation solving, and matrix inverse.
|
inline |
Computes the eigenvalues of the matrix A.
A | a matrix |
evals | on return, the eigenvalues will be stored here in ascending order |
|
inline |
Computes the eigenvalues and eigenvectors of the matrix A.
A | a square symmetric matrix on input, eigenvectors corresponding to eigenvalues on return |
evals | on return, the eigenvalues will be stored here in ascending order |
|
inlinestatic |
Performs the Cholesky factorization of a matrix.
A | the matrix A on input; the factorized (upper triangular) matrix on output |
|
inlinestatic |
Performs the LU factorization of a matrix.
A | the m x n matrix to be factored; on exit, the factors L and U from the factorization M = P*L*U (unit diagonal elements of L are not stored) |
pivwork | on output, contains the pivot indices (for 1 <= i <= min(m,n), row i of A was interchanged with row pivwork[i]) |
|
inline |
Performs the QR factorization of a matrix with column pivoting.
Factorizes A*P = Q*R
AQ | the matrix A on input; the matrix R on output |
Q | the matrix Q on output |
PI | the column pivots on output |
|
inline |
Performs the QR factorization of a matrix.
AQ | the m x n matrix A on input; the matrix min(m,n) x n R on output |
Q | the m x min(m,n) matrix Q on output |
|
static |
Less robust least squares solver (solves Ax = b)
|
inlinestatic |
Inverts the symmetric, positive-definite matrix A using Cholesky factorization.
A | the Cholesky factorization of a matrix; contains the inverse on return |
|
inline |
Computes a matrix inverse using the factorization determined via factor_LU()
M | the LU-factored matrix |
pivwork | the pivots determined by factor_LU() |
|
inlinestatic |
Inverts the symmetric, positive-definite matrix A using Cholesky factorization.
A | a square, symmetric positive-definite matrix; contains the inverse on return |
|
inline |
Inverts the symmetric, indefinite matrix A.
A | a square, symmetric indefinite matrix; inverse will be contained here on return |
|
inline |
Inverts the matrix A using LU factorization.
A | a square matrix; contains the inverse on return |
|
inline |
Computes the nullspace of a matrix.
|
inlinestatic |
Solves a system of linear equations using the factorization determined via factor_chol()
M | the Cholesky decomposition performed by factor_chol() |
XB | the right hand sides on input, the vectors X on return |
|
inline |
Solves the general system AX = B.
A | a square matrix (destroyed on return) |
XB | the matrix B on input, the matrix X on return |
|
inlinestatic |
Solves systems of linear equations using the factorization determined via factor_LDL()
M | the factorization performed by factor_LDL() |
XB | the right hand sides on input, the vectors X on return |
|
inline |
Most robust system of linear equations solver (solves AX = B)
Solves rank-deficient and underdetermined (minimum norm solution) systems. Computes least-squares solution to overdetermined systems.
A | the coefficient matrix (destroyed on return) |
XB | the matrix B on input, the matrix X on return |
tol | the tolerance for determining the rank of A; if tol < 0.0, tol is computed using machine epsilon |
|
inline |
Most robust system of linear equations solver (solves AX = B)
Solves rank-deficient and underdetermined (minimum norm solution) systems. Computes least-squares solution to overdetermined systems.
A | the coefficient matrix (destroyed on return) |
XB | the matrix B on input, the matrix X on return |
tol | the tolerance for determining the rank of A; if tol < 0.0, tol is computed using machine epsilon |
|
inlinestatic |
Solves a system of linear equations using the factorization determined via factor_LU()
M | the LU factorization performed by factor_LU() |
XB | the right hand side on input, the matrix X on return |
transpose | if true, solves M'X = B |
pivwork | pivots computed by factor_LU() |
|
inlinestatic |
Solves a system of equations A*X = B using a symmetric, positive-definite square matrix.
A | the matrix coefficients; this matrix will be destroyed on return |
XB | on input B, on output X |
|
inline |
Solves a symmetric, indefinite square matrix.
A | the matrix to be solved; the matrix is destroyed on return |
XB | the RHS B (A*X = B) on input; the solution, X, on return |
|
inlinestatic |
Solves a triangular system of linear equations.
A | the matrix |
utri | if true A is upper triangular (lower triangular otherwise) |
transpose_A | if true, solves A'*x = b |
XB | contains b on entry, x on return |
|
inlinestatic |
Solves a tridiagonal system.
dl | the (n-1) elements on the subdiagonal (destroyed on return) |
d | the n elements on the diagonal (destroyed on return) |
du | the n elements on the superdiagonal (destroyed on return) |
XB | the right hand side on input, the solution on return |
|
inline |
Does an 'in place' SVD (destroying A)
The singular value decomposition of A is U*S*V' (' is the transpose operator); to recompose A, it will be necessary to transpose V before multiplication (i.e., V is returned by the algorithm, not V'). Note: passed matrices and vectors U, S, and V are resized as necessary.
A | the matrix on which the SVD will be performed (destroyed on return) |
U | on output, a A.rows() x A.rows() orthogonal matrix |
S | on output, a min(A.rows(), A.columns()) length vector of singular values |
V | on output, a A.columns() x A.columns() orthogonal matrix |
|
inline |
Does an 'in place' SVD (destroying A), using divide and conquer algorithm.
The singular value decomposition of A is U*S*V' (' is the transpose operator); to recompose A, it will be necessary to transpose V before multiplication (i.e., V is returned by the algorithm, not V'). Note: passed matrices and vectors U, S, and V are resized as necessary.
A | the matrix on which the SVD will be performed (destroyed on return) |
U | on output, a A.rows() x A.rows() orthogonal matrix |
S | on output, a min(A.rows(), A.columns()) length vector of singular values |
V | on output, a A.columns() x A.columns() orthogonal matrix |
Updates a QR factorization by deleting p columns starting at column idx k.
Q | a m x min(m,n) matrix |
R | a min(m,n) x n matrix |
k | the column index to start deleting at p the number of columns to delete |
Updates a QR factorization by deleting a block of rows.
Q | a m x min(m,n) matrix |
R | a min(m,n) x n matrix |
k | the index to start deleting at |
p | the number of rows to delete |
Updates a QR factorization by inserting one or more columns at column idx k.
Q | a m x min(m,n) matrix |
R | a min(m,n) x n matrix |
U | a m x p matrix, destroyed on return |
Updates a QR factorization by inserting a block of rows, starting at index k.
Q | a m x min(m,n) matrix |
R | a min(m,n) x n matrix |
U | a p x n matrix (destroyed on return) |
k | the index to insert at |
Updates a QR factorization by a rank-1 update.
Q | a m x min(m,n) matrix |
R | a min(m,n) x n matrix |